xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact.c (revision 6df5ee2e43ea7dd5996a6bc97aef5892167e0b91)
149b5e25fSSatish Balay 
2c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
3c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h>
4af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h>
5c6db04a5SJed Brown #include <petscis.h>
649b5e25fSSatish Balay 
75f9f512dSHong Zhang /*
85f9f512dSHong Zhang   input:
9c037c3f7SHong Zhang    F -- numeric factor
105f9f512dSHong Zhang   output:
11c037c3f7SHong Zhang    nneg, nzero, npos: matrix inertia
125f9f512dSHong Zhang */
135f9f512dSHong Zhang 
145f9f512dSHong Zhang #undef __FUNCT__
155f9f512dSHong Zhang #define __FUNCT__ "MatGetInertia_SeqSBAIJ"
1613f74950SBarry Smith PetscErrorCode MatGetInertia_SeqSBAIJ(Mat F,PetscInt *nneig,PetscInt *nzero,PetscInt *npos)
175f9f512dSHong Zhang {
18638f5ce0SDinesh Kaushik   Mat_SeqSBAIJ *fact_ptr = (Mat_SeqSBAIJ*)F->data;
19dd6ea824SBarry Smith   MatScalar    *dd       = fact_ptr->a;
20e70ad169SHong Zhang   PetscInt     mbs       =fact_ptr->mbs,bs=F->rmap->bs,i,nneig_tmp,npos_tmp,*fi = fact_ptr->diag;
215f9f512dSHong Zhang 
225f9f512dSHong Zhang   PetscFunctionBegin;
23e32f2f54SBarry Smith   if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for bs: %D >1 yet",bs);
24eeeff2ecSHong Zhang   nneig_tmp = 0; npos_tmp = 0;
25eeeff2ecSHong Zhang   for (i=0; i<mbs; i++) {
2626fbe8dcSKarl Rupp     if (PetscRealPart(dd[*fi]) > 0.0) npos_tmp++;
2726fbe8dcSKarl Rupp     else if (PetscRealPart(dd[*fi]) < 0.0) nneig_tmp++;
28eeeff2ecSHong Zhang     fi++;
293e0d88b5SBarry Smith   }
30eeeff2ecSHong Zhang   if (nneig) *nneig = nneig_tmp;
31eeeff2ecSHong Zhang   if (npos)  *npos  = npos_tmp;
32eeeff2ecSHong Zhang   if (nzero) *nzero = mbs - nneig_tmp - npos_tmp;
335f9f512dSHong Zhang   PetscFunctionReturn(0);
345f9f512dSHong Zhang }
355f9f512dSHong Zhang 
365f9f512dSHong Zhang /*
375f9f512dSHong Zhang   Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP.
3810c27e3dSHong Zhang   Use Modified Sparse Row (MSR) storage for u and ju. See page 85, "Iterative Methods ..." by Saad.
395f9f512dSHong Zhang */
4010c27e3dSHong Zhang #undef __FUNCT__
4110c27e3dSHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ_MSR"
420481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(Mat F,Mat A,IS perm,const MatFactorInfo *info)
4310c27e3dSHong Zhang {
4410c27e3dSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b;
4510c27e3dSHong Zhang   PetscErrorCode ierr;
465d0c19d7SBarry Smith   const PetscInt *rip,*ai,*aj;
475d0c19d7SBarry Smith   PetscInt       i,mbs = a->mbs,*jutmp,bs = A->rmap->bs,bs2=a->bs2;
4810c27e3dSHong Zhang   PetscInt       m,reallocs = 0,prow;
4910c27e3dSHong Zhang   PetscInt       *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd;
5010c27e3dSHong Zhang   PetscReal      f = info->fill;
51ace3abfcSBarry Smith   PetscBool      perm_identity;
5210c27e3dSHong Zhang 
5310c27e3dSHong Zhang   PetscFunctionBegin;
5410c27e3dSHong Zhang   /* check whether perm is the identity mapping */
5510c27e3dSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
5610c27e3dSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
5710c27e3dSHong Zhang 
5810c27e3dSHong Zhang   if (perm_identity) { /* without permutation */
5910c27e3dSHong Zhang     a->permute = PETSC_FALSE;
6026fbe8dcSKarl Rupp 
6110c27e3dSHong Zhang     ai = a->i; aj = a->j;
6210c27e3dSHong Zhang   } else {            /* non-trivial permutation */
6310c27e3dSHong Zhang     a->permute = PETSC_TRUE;
6426fbe8dcSKarl Rupp 
6510c27e3dSHong Zhang     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
6626fbe8dcSKarl Rupp 
6710c27e3dSHong Zhang     ai = a->inew; aj = a->jnew;
6810c27e3dSHong Zhang   }
6910c27e3dSHong Zhang 
7010c27e3dSHong Zhang   /* initialization */
71854ce69bSBarry Smith   ierr  = PetscMalloc1(mbs+1,&iu);CHKERRQ(ierr);
7210c27e3dSHong Zhang   umax  = (PetscInt)(f*ai[mbs] + 1); umax += mbs + 1;
73785e854fSJed Brown   ierr  = PetscMalloc1(umax,&ju);CHKERRQ(ierr);
7410c27e3dSHong Zhang   iu[0] = mbs+1;
7510c27e3dSHong Zhang   juidx = mbs + 1; /* index for ju */
76d8c74875SBarry Smith   /* jl linked list for pivot row -- linked list for col index */
77dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&jl,mbs,&q);CHKERRQ(ierr);
7810c27e3dSHong Zhang   for (i=0; i<mbs; i++) {
7910c27e3dSHong Zhang     jl[i] = mbs;
8010c27e3dSHong Zhang     q[i]  = 0;
8110c27e3dSHong Zhang   }
8210c27e3dSHong Zhang 
8310c27e3dSHong Zhang   /* for each row k */
8410c27e3dSHong Zhang   for (k=0; k<mbs; k++) {
8510c27e3dSHong Zhang     for (i=0; i<mbs; i++) q[i] = 0;  /* to be removed! */
8610c27e3dSHong Zhang     nzk  = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */
8710c27e3dSHong Zhang     q[k] = mbs;
8810c27e3dSHong Zhang     /* initialize nonzero structure of k-th row to row rip[k] of A */
8910c27e3dSHong Zhang     jmin = ai[rip[k]] +1; /* exclude diag[k] */
9010c27e3dSHong Zhang     jmax = ai[rip[k]+1];
9110c27e3dSHong Zhang     for (j=jmin; j<jmax; j++) {
9210c27e3dSHong Zhang       vj = rip[aj[j]]; /* col. value */
9310c27e3dSHong Zhang       if (vj > k) {
9410c27e3dSHong Zhang         qm = k;
9510c27e3dSHong Zhang         do {
9610c27e3dSHong Zhang           m = qm; qm = q[m];
9710c27e3dSHong Zhang         } while (qm < vj);
98535b19f3SBarry Smith         if (qm == vj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Duplicate entry in A\n");
9910c27e3dSHong Zhang         nzk++;
10010c27e3dSHong Zhang         q[m]  = vj;
10110c27e3dSHong Zhang         q[vj] = qm;
10210c27e3dSHong Zhang       } /* if (vj > k) */
10310c27e3dSHong Zhang     } /* for (j=jmin; j<jmax; j++) */
10410c27e3dSHong Zhang 
10510c27e3dSHong Zhang     /* modify nonzero structure of k-th row by computing fill-in
10610c27e3dSHong Zhang        for each row i to be merged in */
10710c27e3dSHong Zhang     prow = k;
10810c27e3dSHong Zhang     prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */
10910c27e3dSHong Zhang 
11010c27e3dSHong Zhang     while (prow < k) {
11110c27e3dSHong Zhang       /* merge row prow into k-th row */
11210c27e3dSHong Zhang       jmin = iu[prow] + 1; jmax = iu[prow+1];
11310c27e3dSHong Zhang       qm   = k;
11410c27e3dSHong Zhang       for (j=jmin; j<jmax; j++) {
11510c27e3dSHong Zhang         vj = ju[j];
11610c27e3dSHong Zhang         do {
11710c27e3dSHong Zhang           m = qm; qm = q[m];
11810c27e3dSHong Zhang         } while (qm < vj);
11910c27e3dSHong Zhang         if (qm != vj) {
12010c27e3dSHong Zhang           nzk++; q[m] = vj; q[vj] = qm; qm = vj;
12110c27e3dSHong Zhang         }
12210c27e3dSHong Zhang       }
12310c27e3dSHong Zhang       prow = jl[prow]; /* next pivot row */
12410c27e3dSHong Zhang     }
12510c27e3dSHong Zhang 
12610c27e3dSHong Zhang     /* add k to row list for first nonzero element in k-th row */
12710c27e3dSHong Zhang     if (nzk > 0) {
12810c27e3dSHong Zhang       i     = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */
12910c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
13010c27e3dSHong Zhang     }
13110c27e3dSHong Zhang     iu[k+1] = iu[k] + nzk;
13210c27e3dSHong Zhang 
13310c27e3dSHong Zhang     /* allocate more space to ju if needed */
13410c27e3dSHong Zhang     if (iu[k+1] > umax) {
13510c27e3dSHong Zhang       /* estimate how much additional space we will need */
13610c27e3dSHong Zhang       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
13710c27e3dSHong Zhang       /* just double the memory each time */
13810c27e3dSHong Zhang       maxadd = umax;
13910c27e3dSHong Zhang       if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2;
14010c27e3dSHong Zhang       umax += maxadd;
14110c27e3dSHong Zhang 
14210c27e3dSHong Zhang       /* allocate a longer ju */
143785e854fSJed Brown       ierr = PetscMalloc1(umax,&jutmp);CHKERRQ(ierr);
14410c27e3dSHong Zhang       ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(PetscInt));CHKERRQ(ierr);
14510c27e3dSHong Zhang       ierr = PetscFree(ju);CHKERRQ(ierr);
14610c27e3dSHong Zhang       ju   = jutmp;
14710c27e3dSHong Zhang       reallocs++; /* count how many times we realloc */
14810c27e3dSHong Zhang     }
14910c27e3dSHong Zhang 
15010c27e3dSHong Zhang     /* save nonzero structure of k-th row in ju */
15110c27e3dSHong Zhang     i=k;
15210c27e3dSHong Zhang     while (nzk--) {
15310c27e3dSHong Zhang       i           = q[i];
15410c27e3dSHong Zhang       ju[juidx++] = i;
15510c27e3dSHong Zhang     }
15610c27e3dSHong Zhang   }
15710c27e3dSHong Zhang 
1586cf91177SBarry Smith #if defined(PETSC_USE_INFO)
15910c27e3dSHong Zhang   if (ai[mbs] != 0) {
16010c27e3dSHong Zhang     PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
16157622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);CHKERRQ(ierr);
16257622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
16357622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);CHKERRQ(ierr);
164ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
16510c27e3dSHong Zhang   } else {
166ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
16710c27e3dSHong Zhang   }
16863ba0a88SBarry Smith #endif
16910c27e3dSHong Zhang 
17010c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
171d8c74875SBarry Smith   ierr = PetscFree2(jl,q);CHKERRQ(ierr);
17210c27e3dSHong Zhang 
17310c27e3dSHong Zhang   /* put together the new matrix */
1740298fd71SBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(F,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
17510c27e3dSHong Zhang 
1763bb1ff40SBarry Smith   /* ierr = PetscLogObjectParent((PetscObject)B,(PetscObject)iperm);CHKERRQ(ierr); */
177719d5645SBarry Smith   b                = (Mat_SeqSBAIJ*)(F)->data;
17810c27e3dSHong Zhang   b->singlemalloc  = PETSC_FALSE;
179e6b907acSBarry Smith   b->free_a        = PETSC_TRUE;
180e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
18126fbe8dcSKarl Rupp 
182854ce69bSBarry Smith   ierr    = PetscMalloc1((iu[mbs]+1)*bs2,&b->a);CHKERRQ(ierr);
18310c27e3dSHong Zhang   b->j    = ju;
18410c27e3dSHong Zhang   b->i    = iu;
18510c27e3dSHong Zhang   b->diag = 0;
18610c27e3dSHong Zhang   b->ilen = 0;
18710c27e3dSHong Zhang   b->imax = 0;
18810c27e3dSHong Zhang   b->row  = perm;
18926fbe8dcSKarl Rupp 
19010c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
19126fbe8dcSKarl Rupp 
19210c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
19326fbe8dcSKarl Rupp 
19410c27e3dSHong Zhang   b->icol = perm;
19510c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
196854ce69bSBarry Smith   ierr    = PetscMalloc1(bs*mbs+bs,&b->solve_work);CHKERRQ(ierr);
19710c27e3dSHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.
19810c27e3dSHong Zhang      Allocate idnew, solve_work, new a, new j */
1993bb1ff40SBarry Smith   ierr     = PetscLogObjectMemory((PetscObject)F,(iu[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
20010c27e3dSHong Zhang   b->maxnz = b->nz = iu[mbs];
20110c27e3dSHong Zhang 
202719d5645SBarry Smith   (F)->info.factor_mallocs   = reallocs;
203719d5645SBarry Smith   (F)->info.fill_ratio_given = f;
20410c27e3dSHong Zhang   if (ai[mbs] != 0) {
205719d5645SBarry Smith     (F)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
20610c27e3dSHong Zhang   } else {
207719d5645SBarry Smith     (F)->info.fill_ratio_needed = 0.0;
20810c27e3dSHong Zhang   }
209d595f711SHong Zhang   ierr = MatSeqSBAIJSetNumericFactorization_inplace(F,perm_identity);CHKERRQ(ierr);
21010c27e3dSHong Zhang   PetscFunctionReturn(0);
21110c27e3dSHong Zhang }
21210c27e3dSHong Zhang /*
21310c27e3dSHong Zhang     Symbolic U^T*D*U factorization for SBAIJ format.
21480d3d5a7SHong Zhang     See MatICCFactorSymbolic_SeqAIJ() for description of its data structure.
21510c27e3dSHong Zhang */
216c6db04a5SJed Brown #include <petscbt.h>
217c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
2184a2ae208SSatish Balay #undef __FUNCT__
2194a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ"
2200481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
22149b5e25fSSatish Balay {
22298a8e78dSHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data;
22398a8e78dSHong Zhang   Mat_SeqSBAIJ       *b;
2246849ba73SBarry Smith   PetscErrorCode     ierr;
225ace3abfcSBarry Smith   PetscBool          perm_identity,missing;
22698a8e78dSHong Zhang   PetscReal          fill = info->fill;
2277b056e98SHong Zhang   const PetscInt     *rip,*ai=a->i,*aj=a->j;
2289186b105SHong Zhang   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow;
22998a8e78dSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
230c6d0d4f0SHong Zhang   PetscInt           nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr,*udiag;
2310298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
23298a8e78dSHong Zhang   PetscBT            lnkbt;
233d595f711SHong Zhang 
234d595f711SHong Zhang   PetscFunctionBegin;
2359186b105SHong Zhang   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n);
2369186b105SHong Zhang   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
2379186b105SHong Zhang   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
2386d07c2b1SHong Zhang   if (bs > 1) {
2396d07c2b1SHong Zhang     ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(fact,A,perm,info);CHKERRQ(ierr);
2406d07c2b1SHong Zhang     PetscFunctionReturn(0);
2416d07c2b1SHong Zhang   }
242d595f711SHong Zhang 
243d595f711SHong Zhang   /* check whether perm is the identity mapping */
244d595f711SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
245e32f2f54SBarry Smith   if (!perm_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format");
246d595f711SHong Zhang   a->permute = PETSC_FALSE;
247d595f711SHong Zhang   ierr       = ISGetIndices(perm,&rip);CHKERRQ(ierr);
248d595f711SHong Zhang 
249d595f711SHong Zhang   /* initialization */
250854ce69bSBarry Smith   ierr  = PetscMalloc1(mbs+1,&ui);CHKERRQ(ierr);
251854ce69bSBarry Smith   ierr  = PetscMalloc1(mbs+1,&udiag);CHKERRQ(ierr);
252d595f711SHong Zhang   ui[0] = 0;
253d595f711SHong Zhang 
254d595f711SHong Zhang   /* jl: linked list for storing indices of the pivot rows
255d595f711SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
256dcca6d9dSJed Brown   ierr = PetscMalloc4(mbs,&ui_ptr,mbs,&il,mbs,&jl,mbs,&cols);CHKERRQ(ierr);
257d595f711SHong Zhang   for (i=0; i<mbs; i++) {
258d595f711SHong Zhang     jl[i] = mbs; il[i] = 0;
259d595f711SHong Zhang   }
260d595f711SHong Zhang 
261d595f711SHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
262d595f711SHong Zhang   nlnk = mbs + 1;
263d595f711SHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
264d595f711SHong Zhang 
265d595f711SHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
266f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,ai[mbs]+1),&free_space);CHKERRQ(ierr);
267d595f711SHong Zhang   current_space = free_space;
268d595f711SHong Zhang 
269d595f711SHong Zhang   for (k=0; k<mbs; k++) {  /* for each active row k */
270d595f711SHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
271d595f711SHong Zhang     nzk   = 0;
272c6d0d4f0SHong Zhang     ncols = ai[k+1] - ai[k];
273e32f2f54SBarry Smith     if (!ncols) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Empty row %D in matrix ",k);
274d595f711SHong Zhang     for (j=0; j<ncols; j++) {
275c6d0d4f0SHong Zhang       i       = *(aj + ai[k] + j);
276c6d0d4f0SHong Zhang       cols[j] = i;
277d595f711SHong Zhang     }
278d595f711SHong Zhang     ierr = PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
279d595f711SHong Zhang     nzk += nlnk;
280d595f711SHong Zhang 
281d595f711SHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
282d595f711SHong Zhang     prow = jl[k]; /* 1st pivot row */
283d595f711SHong Zhang 
284d595f711SHong Zhang     while (prow < k) {
285d595f711SHong Zhang       nextprow = jl[prow];
286d595f711SHong Zhang       /* merge prow into k-th row */
287d595f711SHong Zhang       jmin   = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
288d595f711SHong Zhang       jmax   = ui[prow+1];
289d595f711SHong Zhang       ncols  = jmax-jmin;
290d595f711SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
291d595f711SHong Zhang       ierr   = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
292d595f711SHong Zhang       nzk   += nlnk;
293d595f711SHong Zhang 
294d595f711SHong Zhang       /* update il and jl for prow */
295d595f711SHong Zhang       if (jmin < jmax) {
296d595f711SHong Zhang         il[prow] = jmin;
297d595f711SHong Zhang         j        = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
298d595f711SHong Zhang       }
299d595f711SHong Zhang       prow = nextprow;
300d595f711SHong Zhang     }
301d595f711SHong Zhang 
302d595f711SHong Zhang     /* if free space is not available, make more free space */
303d595f711SHong Zhang     if (current_space->local_remaining<nzk) {
304d595f711SHong Zhang       i    = mbs - k + 1; /* num of unfactored rows */
305f91af8c7SBarry Smith       i    = PetscIntMultTruncate(i,PetscMin(nzk, i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
306d595f711SHong Zhang       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
307d595f711SHong Zhang       reallocs++;
308d595f711SHong Zhang     }
309d595f711SHong Zhang 
310d595f711SHong Zhang     /* copy data into free space, then initialize lnk */
311d595f711SHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
312d595f711SHong Zhang 
313d595f711SHong Zhang     /* add the k-th row into il and jl */
3147b056e98SHong Zhang     if (nzk > 1) {
315d595f711SHong Zhang       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
316d595f711SHong Zhang       jl[k] = jl[i]; jl[i] = k;
317d595f711SHong Zhang       il[k] = ui[k] + 1;
318d595f711SHong Zhang     }
319d595f711SHong Zhang     ui_ptr[k] = current_space->array;
32026fbe8dcSKarl Rupp 
321d595f711SHong Zhang     current_space->array           += nzk;
322d595f711SHong Zhang     current_space->local_used      += nzk;
323d595f711SHong Zhang     current_space->local_remaining -= nzk;
324d595f711SHong Zhang 
325d595f711SHong Zhang     ui[k+1] = ui[k] + nzk;
326d595f711SHong Zhang   }
327d595f711SHong Zhang 
328d595f711SHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
329d595f711SHong Zhang   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
330d595f711SHong Zhang 
331d595f711SHong Zhang   /* destroy list of free space and other temporary array(s) */
332854ce69bSBarry Smith   ierr = PetscMalloc1(ui[mbs]+1,&uj);CHKERRQ(ierr);
3337b056e98SHong Zhang   ierr = PetscFreeSpaceContiguous_Cholesky(&free_space,uj,mbs,ui,udiag);CHKERRQ(ierr); /* store matrix factor */
334d595f711SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
335d595f711SHong Zhang 
336d595f711SHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
3370298fd71SBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
338d595f711SHong Zhang 
3397b056e98SHong Zhang   b               = (Mat_SeqSBAIJ*)fact->data;
340d595f711SHong Zhang   b->singlemalloc = PETSC_FALSE;
341d595f711SHong Zhang   b->free_a       = PETSC_TRUE;
342d595f711SHong Zhang   b->free_ij      = PETSC_TRUE;
34326fbe8dcSKarl Rupp 
344854ce69bSBarry Smith   ierr = PetscMalloc1(ui[mbs]+1,&b->a);CHKERRQ(ierr);
34526fbe8dcSKarl Rupp 
346d595f711SHong Zhang   b->j         = uj;
347d595f711SHong Zhang   b->i         = ui;
348c6d0d4f0SHong Zhang   b->diag      = udiag;
349c6d0d4f0SHong Zhang   b->free_diag = PETSC_TRUE;
350d595f711SHong Zhang   b->ilen      = 0;
351d595f711SHong Zhang   b->imax      = 0;
352d595f711SHong Zhang   b->row       = perm;
353d595f711SHong Zhang   b->icol      = perm;
35426fbe8dcSKarl Rupp 
355d595f711SHong Zhang   ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
3567b056e98SHong Zhang   ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
35726fbe8dcSKarl Rupp 
3587b056e98SHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
35926fbe8dcSKarl Rupp 
360854ce69bSBarry Smith   ierr = PetscMalloc1(mbs+1,&b->solve_work);CHKERRQ(ierr);
3613bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)fact,ui[mbs]*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
36226fbe8dcSKarl Rupp 
363d595f711SHong Zhang   b->maxnz = b->nz = ui[mbs];
364d595f711SHong Zhang 
36595b5ac22SHong Zhang   fact->info.factor_mallocs   = reallocs;
36695b5ac22SHong Zhang   fact->info.fill_ratio_given = fill;
367d595f711SHong Zhang   if (ai[mbs] != 0) {
36895b5ac22SHong Zhang     fact->info.fill_ratio_needed = ((PetscReal)ui[mbs])/ai[mbs];
369d595f711SHong Zhang   } else {
37095b5ac22SHong Zhang     fact->info.fill_ratio_needed = 0.0;
371d595f711SHong Zhang   }
37295b5ac22SHong Zhang #if defined(PETSC_USE_INFO)
37395b5ac22SHong Zhang   if (ai[mbs] != 0) {
37495b5ac22SHong Zhang     PetscReal af = fact->info.fill_ratio_needed;
37557622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);CHKERRQ(ierr);
37657622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
37757622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);CHKERRQ(ierr);
37895b5ac22SHong Zhang   } else {
37995b5ac22SHong Zhang     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
38095b5ac22SHong Zhang   }
38195b5ac22SHong Zhang #endif
382c6d0d4f0SHong Zhang   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
383d595f711SHong Zhang   PetscFunctionReturn(0);
384d595f711SHong Zhang }
385d595f711SHong Zhang 
386d595f711SHong Zhang #undef __FUNCT__
387d595f711SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ_inplace"
388d595f711SHong Zhang PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
389d595f711SHong Zhang {
390d595f711SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data;
391d595f711SHong Zhang   Mat_SeqSBAIJ       *b;
392d595f711SHong Zhang   PetscErrorCode     ierr;
393ace3abfcSBarry Smith   PetscBool          perm_identity,missing;
394d595f711SHong Zhang   PetscReal          fill = info->fill;
395d595f711SHong Zhang   const PetscInt     *rip,*ai,*aj;
396d595f711SHong Zhang   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow,d;
397d595f711SHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
398d595f711SHong Zhang   PetscInt           nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr;
3990298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
400d595f711SHong Zhang   PetscBT            lnkbt;
40149b5e25fSSatish Balay 
40249b5e25fSSatish Balay   PetscFunctionBegin;
40358ebbce7SBarry Smith   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
404e32f2f54SBarry Smith   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
40558ebbce7SBarry Smith 
40610c27e3dSHong Zhang   /*
40710c27e3dSHong Zhang    This code originally uses Modified Sparse Row (MSR) storage
40810c27e3dSHong Zhang    (see page 85, "Iterative Methods ..." by Saad) for the output matrix B - bad choise!
40910c27e3dSHong Zhang    Then it is rewritten so the factor B takes seqsbaij format. However the associated
41010c27e3dSHong Zhang    MatCholeskyFactorNumeric_() have not been modified for the cases of bs>1 or !perm_identity,
41110c27e3dSHong Zhang    thus the original code in MSR format is still used for these cases.
41210c27e3dSHong Zhang    The code below should replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR() whenever
41310c27e3dSHong Zhang    MatCholeskyFactorNumeric_() is modified for using sbaij symbolic factor.
41410c27e3dSHong Zhang   */
415fff829cfSHong Zhang   if (bs > 1) {
416719d5645SBarry Smith     ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(fact,A,perm,info);CHKERRQ(ierr);
417fff829cfSHong Zhang     PetscFunctionReturn(0);
418fff829cfSHong Zhang   }
41910c27e3dSHong Zhang 
42098a8e78dSHong Zhang   /* check whether perm is the identity mapping */
42198a8e78dSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
42298a8e78dSHong Zhang 
423fff829cfSHong Zhang   if (perm_identity) {
424fff829cfSHong Zhang     a->permute = PETSC_FALSE;
42526fbe8dcSKarl Rupp 
426fff829cfSHong Zhang     ai = a->i; aj = a->j;
4276c4ed002SBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format");
428fff829cfSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
429fff829cfSHong Zhang 
430fff829cfSHong Zhang   /* initialization */
431854ce69bSBarry Smith   ierr  = PetscMalloc1(mbs+1,&ui);CHKERRQ(ierr);
43298a8e78dSHong Zhang   ui[0] = 0;
43398a8e78dSHong Zhang 
43498a8e78dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
4357625dc9aSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
436dcca6d9dSJed Brown   ierr = PetscMalloc4(mbs,&ui_ptr,mbs,&il,mbs,&jl,mbs,&cols);CHKERRQ(ierr);
4377625dc9aSHong Zhang   for (i=0; i<mbs; i++) {
4387625dc9aSHong Zhang     jl[i] = mbs; il[i] = 0;
439fff829cfSHong Zhang   }
440fff829cfSHong Zhang 
44198a8e78dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
4427625dc9aSHong Zhang   nlnk = mbs + 1;
4437625dc9aSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
444fff829cfSHong Zhang 
4457625dc9aSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
446f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,ai[mbs]+1),&free_space);CHKERRQ(ierr);
44798a8e78dSHong Zhang   current_space = free_space;
44898a8e78dSHong Zhang 
4497625dc9aSHong Zhang   for (k=0; k<mbs; k++) {  /* for each active row k */
45098a8e78dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
45198a8e78dSHong Zhang     nzk   = 0;
45298a8e78dSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
4537625dc9aSHong Zhang     for (j=0; j<ncols; j++) {
4547625dc9aSHong Zhang       i       = *(aj + ai[rip[k]] + j);
4557625dc9aSHong Zhang       cols[j] = rip[i];
4567625dc9aSHong Zhang     }
4577625dc9aSHong Zhang     ierr = PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
45898a8e78dSHong Zhang     nzk += nlnk;
45998a8e78dSHong Zhang 
46098a8e78dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
46198a8e78dSHong Zhang     prow = jl[k]; /* 1st pivot row */
462fff829cfSHong Zhang 
463fff829cfSHong Zhang     while (prow < k) {
464fff829cfSHong Zhang       nextprow = jl[prow];
46598a8e78dSHong Zhang       /* merge prow into k-th row */
4667625dc9aSHong Zhang       jmin   = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
46798a8e78dSHong Zhang       jmax   = ui[prow+1];
46898a8e78dSHong Zhang       ncols  = jmax-jmin;
4697625dc9aSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
4705a8e39fbSHong Zhang       ierr   = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
47198a8e78dSHong Zhang       nzk   += nlnk;
472fff829cfSHong Zhang 
47398a8e78dSHong Zhang       /* update il and jl for prow */
474fff829cfSHong Zhang       if (jmin < jmax) {
475fff829cfSHong Zhang         il[prow] = jmin;
47626fbe8dcSKarl Rupp 
4777625dc9aSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
478fff829cfSHong Zhang       }
479fff829cfSHong Zhang       prow = nextprow;
480fff829cfSHong Zhang     }
481fff829cfSHong Zhang 
48298a8e78dSHong Zhang     /* if free space is not available, make more free space */
48398a8e78dSHong Zhang     if (current_space->local_remaining<nzk) {
4847625dc9aSHong Zhang       i    = mbs - k + 1; /* num of unfactored rows */
485f91af8c7SBarry Smith       i    = PetscMin(PetscIntMultTruncate(i,nzk), PetscIntMultTruncate(i,i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
486a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
48798a8e78dSHong Zhang       reallocs++;
48898a8e78dSHong Zhang     }
48998a8e78dSHong Zhang 
49098a8e78dSHong Zhang     /* copy data into free space, then initialize lnk */
4917625dc9aSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
49298a8e78dSHong Zhang 
49398a8e78dSHong Zhang     /* add the k-th row into il and jl */
49498a8e78dSHong Zhang     if (nzk-1 > 0) {
4957625dc9aSHong Zhang       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
496fff829cfSHong Zhang       jl[k] = jl[i]; jl[i] = k;
49798a8e78dSHong Zhang       il[k] = ui[k] + 1;
498fff829cfSHong Zhang     }
4997625dc9aSHong Zhang     ui_ptr[k] = current_space->array;
50026fbe8dcSKarl Rupp 
50198a8e78dSHong Zhang     current_space->array           += nzk;
50298a8e78dSHong Zhang     current_space->local_used      += nzk;
50398a8e78dSHong Zhang     current_space->local_remaining -= nzk;
504fff829cfSHong Zhang 
50598a8e78dSHong Zhang     ui[k+1] = ui[k] + nzk;
506fff829cfSHong Zhang   }
507fff829cfSHong Zhang 
508fff829cfSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
509d8c74875SBarry Smith   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
510fff829cfSHong Zhang 
51198a8e78dSHong Zhang   /* destroy list of free space and other temporary array(s) */
512854ce69bSBarry Smith   ierr = PetscMalloc1(ui[mbs]+1,&uj);CHKERRQ(ierr);
513a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
51498a8e78dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
515fff829cfSHong Zhang 
51698a8e78dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
5170298fd71SBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
51898a8e78dSHong Zhang 
51995b5ac22SHong Zhang   b               = (Mat_SeqSBAIJ*)fact->data;
520fff829cfSHong Zhang   b->singlemalloc = PETSC_FALSE;
521e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
522e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
52326fbe8dcSKarl Rupp 
524854ce69bSBarry Smith   ierr = PetscMalloc1(ui[mbs]+1,&b->a);CHKERRQ(ierr);
52526fbe8dcSKarl Rupp 
52698a8e78dSHong Zhang   b->j    = uj;
52798a8e78dSHong Zhang   b->i    = ui;
528fff829cfSHong Zhang   b->diag = 0;
529fff829cfSHong Zhang   b->ilen = 0;
530fff829cfSHong Zhang   b->imax = 0;
531fff829cfSHong Zhang   b->row  = perm;
53226fbe8dcSKarl Rupp 
533fff829cfSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
53426fbe8dcSKarl Rupp 
535fff829cfSHong Zhang   ierr     = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
536fff829cfSHong Zhang   b->icol  = perm;
537fff829cfSHong Zhang   ierr     = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
538854ce69bSBarry Smith   ierr     = PetscMalloc1(mbs+1,&b->solve_work);CHKERRQ(ierr);
5393bb1ff40SBarry Smith   ierr     = PetscLogObjectMemory((PetscObject)fact,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
5407625dc9aSHong Zhang   b->maxnz = b->nz = ui[mbs];
541fff829cfSHong Zhang 
54295b5ac22SHong Zhang   fact->info.factor_mallocs   = reallocs;
54395b5ac22SHong Zhang   fact->info.fill_ratio_given = fill;
5447625dc9aSHong Zhang   if (ai[mbs] != 0) {
54595b5ac22SHong Zhang     fact->info.fill_ratio_needed = ((PetscReal)ui[mbs])/ai[mbs];
546fff829cfSHong Zhang   } else {
54795b5ac22SHong Zhang     fact->info.fill_ratio_needed = 0.0;
548fff829cfSHong Zhang   }
54995b5ac22SHong Zhang #if defined(PETSC_USE_INFO)
55095b5ac22SHong Zhang   if (ai[mbs] != 0) {
55195b5ac22SHong Zhang     PetscReal af = fact->info.fill_ratio_needed;
55257622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);CHKERRQ(ierr);
55357622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
55457622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);CHKERRQ(ierr);
55595b5ac22SHong Zhang   } else {
55695b5ac22SHong Zhang     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
55795b5ac22SHong Zhang   }
55895b5ac22SHong Zhang #endif
559d595f711SHong Zhang   ierr = MatSeqSBAIJSetNumericFactorization_inplace(fact,perm_identity);CHKERRQ(ierr);
560064503c5SHong Zhang   PetscFunctionReturn(0);
561064503c5SHong Zhang }
562d595f711SHong Zhang 
5634a2ae208SSatish Balay #undef __FUNCT__
5644a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N"
5650481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
56649b5e25fSSatish Balay {
5674c16a6a6SHong Zhang   Mat_SeqSBAIJ   *a   = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
5684c16a6a6SHong Zhang   IS             perm = b->row;
5696849ba73SBarry Smith   PetscErrorCode ierr;
5705d0c19d7SBarry Smith   const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j;
5715d0c19d7SBarry Smith   PetscInt       i,j;
5725d0c19d7SBarry Smith   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
5733bc0b13bSBarry Smith   PetscInt       bs  =A->rmap->bs,bs2 = a->bs2;
5744c16a6a6SHong Zhang   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
5754c16a6a6SHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
57628de702eSHong Zhang   MatScalar      *work;
57713f74950SBarry Smith   PetscInt       *pivots;
5785f8bbccaSHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
5794c16a6a6SHong Zhang 
5804c16a6a6SHong Zhang   PetscFunctionBegin;
5814c16a6a6SHong Zhang   /* initialization */
5821795a4d1SJed Brown   ierr = PetscCalloc1(bs2*mbs,&rtmp);CHKERRQ(ierr);
583dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
5845f8bbccaSHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
5855f8bbccaSHong Zhang 
586*6df5ee2eSHong Zhang   il[0] = 0;
587*6df5ee2eSHong Zhang   for (i=0; i<mbs; i++) jl[i] = mbs;
588*6df5ee2eSHong Zhang 
589dcca6d9dSJed Brown   ierr = PetscMalloc3(bs2,&dk,bs2,&uik,bs,&work);CHKERRQ(ierr);
590785e854fSJed Brown   ierr = PetscMalloc1(bs,&pivots);CHKERRQ(ierr);
5914c16a6a6SHong Zhang 
5924c16a6a6SHong Zhang   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
5934c16a6a6SHong Zhang 
5944c16a6a6SHong Zhang   /* check permutation */
5954c16a6a6SHong Zhang   if (!a->permute) {
5964c16a6a6SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
5974c16a6a6SHong Zhang   } else {
5984c16a6a6SHong Zhang     ai   = a->inew; aj = a->jnew;
599785e854fSJed Brown     ierr = PetscMalloc1(bs2*ai[mbs],&aa);CHKERRQ(ierr);
6004c16a6a6SHong Zhang     ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
601785e854fSJed Brown     ierr = PetscMalloc1(ai[mbs],&a2anew);CHKERRQ(ierr);
60213f74950SBarry Smith     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
6034c16a6a6SHong 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);
65696b95a6bSBarry Smith       PetscKernel_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) */
65996b95a6bSBarry Smith       PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
6603bc0b13bSBarry Smith       ierr = PetscLogFlops(4.0*bs*bs2);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;
67296b95a6bSBarry Smith           PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
6734c16a6a6SHong Zhang         }
6743bc0b13bSBarry Smith         ierr = PetscLogFlops(2.0*bs*bs2*(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);
6895f8bbccaSHong Zhang 
6905f8bbccaSHong Zhang     ierr = PetscKernel_A_gets_inverse_A(bs,diag,pivots,work,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
6915f8bbccaSHong Zhang     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
6924c16a6a6SHong Zhang 
6934c16a6a6SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
6944c16a6a6SHong Zhang     if (jmin < jmax) {
6954c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++) {
6964c16a6a6SHong Zhang         vj       = bj[j];      /* block col. index of U */
6974c16a6a6SHong Zhang         u        = ba + j*bs2;
6984c16a6a6SHong Zhang         rtmp_ptr = rtmp + vj*bs2;
6994c16a6a6SHong Zhang         for (k1=0; k1<bs2; k1++) {
7004c16a6a6SHong Zhang           *u++        = *rtmp_ptr;
7014c16a6a6SHong Zhang           *rtmp_ptr++ = 0.0;
7024c16a6a6SHong Zhang         }
7034c16a6a6SHong Zhang       }
7044c16a6a6SHong Zhang 
7054c16a6a6SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
7064c16a6a6SHong Zhang       il[k] = jmin;
7074c16a6a6SHong Zhang       i     = bj[jmin];
7084c16a6a6SHong Zhang       jl[k] = jl[i]; jl[i] = k;
7094c16a6a6SHong Zhang     }
7104c16a6a6SHong Zhang   }
7114c16a6a6SHong Zhang 
7124c16a6a6SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
713d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
714d8c74875SBarry Smith   ierr = PetscFree3(dk,uik,work);CHKERRQ(ierr);
7154c16a6a6SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
7164c16a6a6SHong Zhang   if (a->permute) {
7174c16a6a6SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
7184c16a6a6SHong Zhang   }
7194c16a6a6SHong Zhang 
7204c16a6a6SHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
72126fbe8dcSKarl Rupp 
7224f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_N_inplace;
7234f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_inplace;
7244f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_inplace;
7254f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_inplace;
726db4efbfdSBarry Smith 
7274c16a6a6SHong Zhang   C->assembled    = PETSC_TRUE;
7284c16a6a6SHong Zhang   C->preallocated = PETSC_TRUE;
72926fbe8dcSKarl Rupp 
730efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
7314c16a6a6SHong Zhang   PetscFunctionReturn(0);
7324c16a6a6SHong Zhang }
733d4edadd4SHong Zhang 
7344a2ae208SSatish Balay #undef __FUNCT__
7354a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering"
7360481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
737671cb588SHong Zhang {
738671cb588SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
739dfbe8321SBarry Smith   PetscErrorCode ierr;
74013f74950SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
74113f74950SBarry Smith   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
7423bc0b13bSBarry Smith   PetscInt       bs  =A->rmap->bs,bs2 = a->bs2;
743671cb588SHong Zhang   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
744671cb588SHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
74528de702eSHong Zhang   MatScalar      *work;
74613f74950SBarry Smith   PetscInt       *pivots;
7475f8bbccaSHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
748671cb588SHong Zhang 
749671cb588SHong Zhang   PetscFunctionBegin;
7501795a4d1SJed Brown   ierr = PetscCalloc1(bs2*mbs,&rtmp);CHKERRQ(ierr);
751dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
752*6df5ee2eSHong Zhang   il[0] = 0;
753*6df5ee2eSHong Zhang   for (i=0; i<mbs; i++) jl[i] = mbs;
754*6df5ee2eSHong Zhang 
755dcca6d9dSJed Brown   ierr = PetscMalloc3(bs2,&dk,bs2,&uik,bs,&work);CHKERRQ(ierr);
756785e854fSJed Brown   ierr = PetscMalloc1(bs,&pivots);CHKERRQ(ierr);
7575f8bbccaSHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
758671cb588SHong Zhang 
759671cb588SHong Zhang   ai = a->i; aj = a->j; aa = a->a;
760671cb588SHong Zhang 
761671cb588SHong Zhang   /* for each row k */
762671cb588SHong Zhang   for (k = 0; k<mbs; k++) {
763671cb588SHong Zhang 
764671cb588SHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
765671cb588SHong Zhang     jmin = ai[k]; jmax = ai[k+1];
766671cb588SHong Zhang     ap   = aa + jmin*bs2;
767671cb588SHong Zhang     for (j = jmin; j < jmax; j++) {
768671cb588SHong Zhang       vj       = aj[j];   /* block col. index */
769671cb588SHong Zhang       rtmp_ptr = rtmp + vj*bs2;
770671cb588SHong Zhang       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
771671cb588SHong Zhang     }
772671cb588SHong Zhang 
773671cb588SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
774671cb588SHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
775671cb588SHong Zhang     i    = jl[k]; /* first row to be added to k_th row  */
776671cb588SHong Zhang 
777057f5ba7SHong Zhang     while (i < k) {
778671cb588SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
779671cb588SHong Zhang 
780671cb588SHong Zhang       /* compute multiplier */
781671cb588SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
782671cb588SHong Zhang 
783671cb588SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
784671cb588SHong Zhang       diag = ba + i*bs2;
785671cb588SHong Zhang       u    = ba + ili*bs2;
786671cb588SHong Zhang       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
78796b95a6bSBarry Smith       PetscKernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
788671cb588SHong Zhang 
789671cb588SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
79096b95a6bSBarry Smith       PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
7913bc0b13bSBarry Smith       ierr = PetscLogFlops(2.0*bs*bs2);CHKERRQ(ierr);
792671cb588SHong Zhang 
793671cb588SHong Zhang       /* update -U(i,k) */
794671cb588SHong Zhang       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
795671cb588SHong Zhang 
796671cb588SHong Zhang       /* add multiple of row i to k-th row ... */
797671cb588SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
798671cb588SHong Zhang       if (jmin < jmax) {
799671cb588SHong Zhang         for (j=jmin; j<jmax; j++) {
800671cb588SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
801671cb588SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
802671cb588SHong Zhang           u        = ba + j*bs2;
80396b95a6bSBarry Smith           PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
804671cb588SHong Zhang         }
8053bc0b13bSBarry Smith         ierr = PetscLogFlops(2.0*bs*bs2*(jmax-jmin));CHKERRQ(ierr);
806671cb588SHong Zhang 
807671cb588SHong Zhang         /* ... add i to row list for next nonzero entry */
808671cb588SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
809671cb588SHong Zhang         j     = bj[jmin];
810671cb588SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
811671cb588SHong Zhang       }
812671cb588SHong Zhang       i = nexti;
813671cb588SHong Zhang     }
814671cb588SHong Zhang 
815671cb588SHong Zhang     /* save nonzero entries in k-th row of U ... */
816671cb588SHong Zhang 
817671cb588SHong Zhang     /* invert diagonal block */
818671cb588SHong Zhang     diag = ba+k*bs2;
819671cb588SHong Zhang     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
8205f8bbccaSHong Zhang 
8215f8bbccaSHong Zhang     ierr = PetscKernel_A_gets_inverse_A(bs,diag,pivots,work,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
8225f8bbccaSHong Zhang     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
823671cb588SHong Zhang 
824671cb588SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
825671cb588SHong Zhang     if (jmin < jmax) {
826671cb588SHong Zhang       for (j=jmin; j<jmax; j++) {
827671cb588SHong Zhang         vj       = bj[j];      /* block col. index of U */
828671cb588SHong Zhang         u        = ba + j*bs2;
829671cb588SHong Zhang         rtmp_ptr = rtmp + vj*bs2;
830671cb588SHong Zhang         for (k1=0; k1<bs2; k1++) {
831671cb588SHong Zhang           *u++        = *rtmp_ptr;
832671cb588SHong Zhang           *rtmp_ptr++ = 0.0;
833671cb588SHong Zhang         }
834671cb588SHong Zhang       }
835671cb588SHong Zhang 
836671cb588SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
837671cb588SHong Zhang       il[k] = jmin;
838671cb588SHong Zhang       i     = bj[jmin];
839671cb588SHong Zhang       jl[k] = jl[i]; jl[i] = k;
840671cb588SHong Zhang     }
841671cb588SHong Zhang   }
842671cb588SHong Zhang 
843671cb588SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
844d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
845d8c74875SBarry Smith   ierr = PetscFree3(dk,uik,work);CHKERRQ(ierr);
846671cb588SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
847671cb588SHong Zhang 
8484f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8494f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8504f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8514f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
852671cb588SHong Zhang   C->assembled           = PETSC_TRUE;
853671cb588SHong Zhang   C->preallocated        = PETSC_TRUE;
85426fbe8dcSKarl Rupp 
855efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
856671cb588SHong Zhang   PetscFunctionReturn(0);
857671cb588SHong Zhang }
858671cb588SHong Zhang 
85949b5e25fSSatish Balay /*
860fcf159c0SHong Zhang     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
861cc0c071aSHong Zhang     Version for blocks 2 by 2.
86249b5e25fSSatish Balay */
8634a2ae208SSatish Balay #undef __FUNCT__
8644a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2"
8650481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat C,Mat A,const MatFactorInfo *info)
86649b5e25fSSatish Balay {
86791602c66SHong Zhang   Mat_SeqSBAIJ   *a   = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
868cc0c071aSHong Zhang   IS             perm = b->row;
8696849ba73SBarry Smith   PetscErrorCode ierr;
8705d0c19d7SBarry Smith   const PetscInt *ai,*aj,*perm_ptr;
8715d0c19d7SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
8725d0c19d7SBarry Smith   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
873d8c74875SBarry Smith   MatScalar      *ba = b->a,*aa,*ap;
874d8c74875SBarry Smith   MatScalar      *u,*diag,*rtmp,*rtmp_ptr,dk[4],uik[4];
87514d2772eSHong Zhang   PetscReal      shift = info->shiftamount;
876a455e926SHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
87749b5e25fSSatish Balay 
87849b5e25fSSatish Balay   PetscFunctionBegin;
8790164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
8800164db54SHong Zhang 
88191602c66SHong Zhang   /* initialization */
88291602c66SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
88391602c66SHong Zhang      window U(0:k, k:mbs-1).
88491602c66SHong Zhang      jl:    list of rows to be added to uneliminated rows
88591602c66SHong Zhang             i>= k: jl(i) is the first row to be added to row i
88691602c66SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
88791602c66SHong Zhang             jl(i) = mbs indicates the end of a list
88891602c66SHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
88991602c66SHong Zhang             row i of U */
8901795a4d1SJed Brown   ierr = PetscCalloc1(4*mbs,&rtmp);CHKERRQ(ierr);
891dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
892*6df5ee2eSHong Zhang   il[0] = 0;
893*6df5ee2eSHong Zhang   for (i=0; i<mbs; i++) jl[i] = mbs;
894*6df5ee2eSHong Zhang 
895cc0c071aSHong Zhang   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
896cc0c071aSHong Zhang 
897cc0c071aSHong Zhang   /* check permutation */
898cc0c071aSHong Zhang   if (!a->permute) {
899cc0c071aSHong Zhang     ai = a->i; aj = a->j; aa = a->a;
900cc0c071aSHong Zhang   } else {
901cc0c071aSHong Zhang     ai   = a->inew; aj = a->jnew;
902785e854fSJed Brown     ierr = PetscMalloc1(4*ai[mbs],&aa);CHKERRQ(ierr);
903cc0c071aSHong Zhang     ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
904785e854fSJed Brown     ierr = PetscMalloc1(ai[mbs],&a2anew);CHKERRQ(ierr);
90513f74950SBarry Smith     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
906cc0c071aSHong Zhang 
907cc0c071aSHong Zhang     for (i=0; i<mbs; i++) {
908cc0c071aSHong Zhang       jmin = ai[i]; jmax = ai[i+1];
909cc0c071aSHong Zhang       for (j=jmin; j<jmax; j++) {
910cc0c071aSHong Zhang         while (a2anew[j] != j) {
911cc0c071aSHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
912cc0c071aSHong Zhang           for (k1=0; k1<4; k1++) {
913cc0c071aSHong Zhang             dk[k1]     = aa[k*4+k1];
914cc0c071aSHong Zhang             aa[k*4+k1] = aa[j*4+k1];
915cc0c071aSHong Zhang             aa[j*4+k1] = dk[k1];
916cc0c071aSHong Zhang           }
917cc0c071aSHong Zhang         }
918cc0c071aSHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
919cc0c071aSHong Zhang         if (i > aj[j]) {
920a1723e09SHong Zhang           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
921cc0c071aSHong Zhang           ap    = aa + j*4;  /* ptr to the beginning of the block */
922cc0c071aSHong Zhang           dk[1] = ap[1];     /* swap ap[1] and ap[2] */
923cc0c071aSHong Zhang           ap[1] = ap[2];
924cc0c071aSHong Zhang           ap[2] = dk[1];
925cc0c071aSHong Zhang         }
926cc0c071aSHong Zhang       }
927cc0c071aSHong Zhang     }
928ac355199SBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
929cc0c071aSHong Zhang   }
9303845f261SHong Zhang 
93191602c66SHong Zhang   /* for each row k */
93291602c66SHong Zhang   for (k = 0; k<mbs; k++) {
93391602c66SHong Zhang 
93491602c66SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
935cc0c071aSHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
936cc0c071aSHong Zhang     ap   = aa + jmin*4;
93791602c66SHong Zhang     for (j = jmin; j < jmax; j++) {
938cc0c071aSHong Zhang       vj       = perm_ptr[aj[j]];   /* block col. index */
939cc0c071aSHong Zhang       rtmp_ptr = rtmp + vj*4;
940cc0c071aSHong Zhang       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
94191602c66SHong Zhang     }
94291602c66SHong Zhang 
94391602c66SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
944cc0c071aSHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
94591602c66SHong Zhang     i    = jl[k]; /* first row to be added to k_th row  */
94691602c66SHong Zhang 
947057f5ba7SHong Zhang     while (i < k) {
94891602c66SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
94991602c66SHong Zhang 
9503845f261SHong Zhang       /* compute multiplier */
95191602c66SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
9523845f261SHong Zhang 
9533845f261SHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
954cc0c071aSHong Zhang       diag   = ba + i*4;
955cc0c071aSHong Zhang       u      = ba + ili*4;
956cc0c071aSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
957cc0c071aSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
958cc0c071aSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
959cc0c071aSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
9603845f261SHong Zhang 
9613845f261SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
962cc0c071aSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
963cc0c071aSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
964cc0c071aSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
965cc0c071aSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
9663845f261SHong Zhang 
967dc0b31edSSatish Balay       ierr = PetscLogFlops(16.0*2.0);CHKERRQ(ierr);
968187a9f4bSHong Zhang 
9693845f261SHong Zhang       /* update -U(i,k): ba[ili] = uik */
970cc0c071aSHong Zhang       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
97191602c66SHong Zhang 
97291602c66SHong Zhang       /* add multiple of row i to k-th row ... */
97391602c66SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
97491602c66SHong Zhang       if (jmin < jmax) {
9753845f261SHong Zhang         for (j=jmin; j<jmax; j++) {
9763845f261SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
977cc0c071aSHong Zhang           rtmp_ptr     = rtmp + bj[j]*4;
978cc0c071aSHong Zhang           u            = ba + j*4;
979cc0c071aSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
980cc0c071aSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
981cc0c071aSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
982cc0c071aSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
9833845f261SHong Zhang         }
984dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*(jmax-jmin));CHKERRQ(ierr);
9853845f261SHong Zhang 
98691602c66SHong Zhang         /* ... add i to row list for next nonzero entry */
98791602c66SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
98891602c66SHong Zhang         j     = bj[jmin];
98991602c66SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
99091602c66SHong Zhang       }
991a1723e09SHong Zhang       i = nexti;
99291602c66SHong Zhang     }
993cc0c071aSHong Zhang 
99491602c66SHong Zhang     /* save nonzero entries in k-th row of U ... */
9953845f261SHong Zhang 
996cc0c071aSHong Zhang     /* invert diagonal block */
997cc0c071aSHong Zhang     diag = ba+k*4;
998cc0c071aSHong Zhang     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
999a455e926SHong Zhang     ierr = PetscKernel_A_gets_inverse_A_2(diag,shift,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
10002e92ee13SHong Zhang     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
10013845f261SHong Zhang 
100291602c66SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
100391602c66SHong Zhang     if (jmin < jmax) {
100491602c66SHong Zhang       for (j=jmin; j<jmax; j++) {
1005cc0c071aSHong Zhang         vj       = bj[j];      /* block col. index of U */
1006cc0c071aSHong Zhang         u        = ba + j*4;
1007cc0c071aSHong Zhang         rtmp_ptr = rtmp + vj*4;
1008cc0c071aSHong Zhang         for (k1=0; k1<4; k1++) {
1009cc0c071aSHong Zhang           *u++        = *rtmp_ptr;
1010cc0c071aSHong Zhang           *rtmp_ptr++ = 0.0;
1011cc0c071aSHong Zhang         }
1012cc0c071aSHong Zhang       }
10133845f261SHong Zhang 
101491602c66SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
101591602c66SHong Zhang       il[k] = jmin;
101691602c66SHong Zhang       i     = bj[jmin];
101791602c66SHong Zhang       jl[k] = jl[i]; jl[i] = k;
101891602c66SHong Zhang     }
101991602c66SHong Zhang   }
10203845f261SHong Zhang 
102149b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1022d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
102391602c66SHong Zhang   if (a->permute) {
102491602c66SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
102591602c66SHong Zhang   }
1026cc0c071aSHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
102726fbe8dcSKarl Rupp 
10284f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_2_inplace;
10294f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_inplace;
103049b5e25fSSatish Balay   C->assembled           = PETSC_TRUE;
10315c0bcdfcSHong Zhang   C->preallocated        = PETSC_TRUE;
103226fbe8dcSKarl Rupp 
1033efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
103449b5e25fSSatish Balay   PetscFunctionReturn(0);
103549b5e25fSSatish Balay }
103691602c66SHong Zhang 
103749b5e25fSSatish Balay /*
103849b5e25fSSatish Balay       Version for when blocks are 2 by 2 Using natural ordering
103949b5e25fSSatish Balay */
10404a2ae208SSatish Balay #undef __FUNCT__
10414a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering"
10420481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
104349b5e25fSSatish Balay {
1044ab27746eSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
1045dfbe8321SBarry Smith   PetscErrorCode ierr;
104613f74950SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
104713f74950SBarry Smith   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
1048d8c74875SBarry Smith   MatScalar      *ba = b->a,*aa,*ap,dk[8],uik[8];
1049ab27746eSHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
105014d2772eSHong Zhang   PetscReal      shift = info->shiftamount;
1051a455e926SHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
105249b5e25fSSatish Balay 
105349b5e25fSSatish Balay   PetscFunctionBegin;
10540164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
10550164db54SHong Zhang 
1056ab27746eSHong Zhang   /* initialization */
1057ab27746eSHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
1058ab27746eSHong Zhang      window U(0:k, k:mbs-1).
1059ab27746eSHong Zhang      jl:    list of rows to be added to uneliminated rows
1060ab27746eSHong Zhang             i>= k: jl(i) is the first row to be added to row i
1061ab27746eSHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
1062ab27746eSHong Zhang             jl(i) = mbs indicates the end of a list
1063ab27746eSHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
1064ab27746eSHong Zhang             row i of U */
10651795a4d1SJed Brown   ierr = PetscCalloc1(4*mbs,&rtmp);CHKERRQ(ierr);
1066dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
1067*6df5ee2eSHong Zhang   il[0] = 0;
1068*6df5ee2eSHong Zhang   for (i=0; i<mbs; i++) jl[i] = mbs;
1069*6df5ee2eSHong Zhang 
1070ab27746eSHong Zhang   ai = a->i; aj = a->j; aa = a->a;
1071ab27746eSHong Zhang 
1072ab27746eSHong Zhang   /* for each row k */
1073ab27746eSHong Zhang   for (k = 0; k<mbs; k++) {
1074ab27746eSHong Zhang 
1075ab27746eSHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
1076ab27746eSHong Zhang     jmin = ai[k]; jmax = ai[k+1];
1077ab27746eSHong Zhang     ap   = aa + jmin*4;
1078ab27746eSHong Zhang     for (j = jmin; j < jmax; j++) {
1079ab27746eSHong Zhang       vj       = aj[j];   /* block col. index */
1080ab27746eSHong Zhang       rtmp_ptr = rtmp + vj*4;
1081ab27746eSHong Zhang       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
108249b5e25fSSatish Balay     }
1083ab27746eSHong Zhang 
1084ab27746eSHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
1085ab27746eSHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
1086ab27746eSHong Zhang     i    = jl[k]; /* first row to be added to k_th row  */
1087ab27746eSHong Zhang 
1088057f5ba7SHong Zhang     while (i < k) {
1089ab27746eSHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
1090ab27746eSHong Zhang 
1091ab27746eSHong Zhang       /* compute multiplier */
1092ab27746eSHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1093ab27746eSHong Zhang 
1094ab27746eSHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
1095ab27746eSHong Zhang       diag   = ba + i*4;
1096ab27746eSHong Zhang       u      = ba + ili*4;
1097ab27746eSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
1098ab27746eSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
1099ab27746eSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
1100ab27746eSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
1101ab27746eSHong Zhang 
1102ab27746eSHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
1103ab27746eSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
1104ab27746eSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
1105ab27746eSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
1106ab27746eSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
1107ab27746eSHong Zhang 
1108dc0b31edSSatish Balay       ierr = PetscLogFlops(16.0*2.0);CHKERRQ(ierr);
1109187a9f4bSHong Zhang 
1110ab27746eSHong Zhang       /* update -U(i,k): ba[ili] = uik */
1111ab27746eSHong Zhang       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
1112ab27746eSHong Zhang 
1113ab27746eSHong Zhang       /* add multiple of row i to k-th row ... */
1114ab27746eSHong Zhang       jmin = ili + 1; jmax = bi[i+1];
1115ab27746eSHong Zhang       if (jmin < jmax) {
1116ab27746eSHong Zhang         for (j=jmin; j<jmax; j++) {
1117ab27746eSHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
1118ab27746eSHong Zhang           rtmp_ptr     = rtmp + bj[j]*4;
1119ab27746eSHong Zhang           u            = ba + j*4;
1120ab27746eSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
1121ab27746eSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
1122ab27746eSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
1123ab27746eSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
112449b5e25fSSatish Balay         }
1125dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*(jmax-jmin));CHKERRQ(ierr);
1126ab27746eSHong Zhang 
1127ab27746eSHong Zhang         /* ... add i to row list for next nonzero entry */
1128ab27746eSHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
1129ab27746eSHong Zhang         j     = bj[jmin];
1130ab27746eSHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
113149b5e25fSSatish Balay       }
1132ab27746eSHong Zhang       i = nexti;
113349b5e25fSSatish Balay     }
1134ab27746eSHong Zhang 
1135ab27746eSHong Zhang     /* save nonzero entries in k-th row of U ... */
1136ab27746eSHong Zhang 
113749b5e25fSSatish Balay     /* invert diagonal block */
1138ab27746eSHong Zhang     diag = ba+k*4;
1139ab27746eSHong Zhang     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
1140a455e926SHong Zhang     ierr = PetscKernel_A_gets_inverse_A_2(diag,shift,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
11412e92ee13SHong Zhang     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1142ab27746eSHong Zhang 
1143ab27746eSHong Zhang     jmin = bi[k]; jmax = bi[k+1];
1144ab27746eSHong Zhang     if (jmin < jmax) {
1145ab27746eSHong Zhang       for (j=jmin; j<jmax; j++) {
1146ab27746eSHong Zhang         vj       = bj[j];      /* block col. index of U */
1147ab27746eSHong Zhang         u        = ba + j*4;
1148ab27746eSHong Zhang         rtmp_ptr = rtmp + vj*4;
1149ab27746eSHong Zhang         for (k1=0; k1<4; k1++) {
1150ab27746eSHong Zhang           *u++        = *rtmp_ptr;
1151ab27746eSHong Zhang           *rtmp_ptr++ = 0.0;
1152ab27746eSHong Zhang         }
1153ab27746eSHong Zhang       }
1154ab27746eSHong Zhang 
1155ab27746eSHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
1156ab27746eSHong Zhang       il[k] = jmin;
1157ab27746eSHong Zhang       i     = bj[jmin];
1158ab27746eSHong Zhang       jl[k] = jl[i]; jl[i] = k;
1159ab27746eSHong Zhang     }
116049b5e25fSSatish Balay   }
116149b5e25fSSatish Balay 
116249b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1163d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
1164ab27746eSHong Zhang 
11654f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11664f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11674f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11684f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
116949b5e25fSSatish Balay   C->assembled           = PETSC_TRUE;
11705c0bcdfcSHong Zhang   C->preallocated        = PETSC_TRUE;
117126fbe8dcSKarl Rupp 
1172efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
117349b5e25fSSatish Balay   PetscFunctionReturn(0);
117449b5e25fSSatish Balay }
117549b5e25fSSatish Balay 
117649b5e25fSSatish Balay /*
117798a8e78dSHong Zhang     Numeric U^T*D*U factorization for SBAIJ format.
117891602c66SHong Zhang     Version for blocks are 1 by 1.
117949b5e25fSSatish Balay */
11804a2ae208SSatish Balay #undef __FUNCT__
1181d595f711SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace"
1182d595f711SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info)
118349b5e25fSSatish Balay {
118449b5e25fSSatish Balay   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ*)C->data;
118549b5e25fSSatish Balay   IS             ip=b->row;
11866849ba73SBarry Smith   PetscErrorCode ierr;
11875d0c19d7SBarry Smith   const PetscInt *ai,*aj,*rip;
11885d0c19d7SBarry Smith   PetscInt       *a2anew,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol;
1189997a0949SHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1190997a0949SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa,dk,uikdi;
11910e95ead3SHong Zhang   PetscReal      rs;
11920e95ead3SHong Zhang   FactorShiftCtx sctx;
119349b5e25fSSatish Balay 
119449b5e25fSSatish Balay   PetscFunctionBegin;
11950e95ead3SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
11960e95ead3SHong Zhang   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
11973cea4cbeSHong Zhang 
119849b5e25fSSatish Balay   ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr);
1199cb718733SHong Zhang   if (!a->permute) {
12002d007305SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
12012d007305SHong Zhang   } else {
12022d007305SHong Zhang     ai     = a->inew; aj = a->jnew;
1203fff829cfSHong Zhang     nz     = ai[mbs];
1204785e854fSJed Brown     ierr   = PetscMalloc1(nz,&aa);CHKERRQ(ierr);
1205fff829cfSHong Zhang     a2anew = a->a2anew;
1206997a0949SHong Zhang     bval   = a->a;
1207fff829cfSHong Zhang     for (j=0; j<nz; j++) {
1208997a0949SHong Zhang       aa[a2anew[j]] = *(bval++);
12092d007305SHong Zhang     }
12102d007305SHong Zhang   }
12112d007305SHong Zhang 
121291602c66SHong Zhang   /* initialization */
121349b5e25fSSatish Balay   /* il and jl record the first nonzero element in each row of the accessing
121449b5e25fSSatish Balay      window U(0:k, k:mbs-1).
121549b5e25fSSatish Balay      jl:    list of rows to be added to uneliminated rows
121649b5e25fSSatish Balay             i>= k: jl(i) is the first row to be added to row i
121749b5e25fSSatish Balay             i<  k: jl(i) is the row following row i in some list of rows
121849b5e25fSSatish Balay             jl(i) = mbs indicates the end of a list
121949b5e25fSSatish Balay      il(i): points to the first nonzero element in columns k,...,mbs-1 of
122049b5e25fSSatish Balay             row i of U */
1221dcca6d9dSJed Brown   ierr = PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&jl);CHKERRQ(ierr);
1222997a0949SHong Zhang 
1223997a0949SHong Zhang   do {
122407b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
1225*6df5ee2eSHong Zhang     il[0] = 0;
122649b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
1227*6df5ee2eSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs;
122849b5e25fSSatish Balay     }
1229997a0949SHong Zhang 
123049b5e25fSSatish Balay     for (k = 0; k<mbs; k++) {
1231997a0949SHong Zhang       /*initialize k-th row by the perm[k]-th row of A */
123249b5e25fSSatish Balay       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
12337701de02SHong Zhang       bval = ba + bi[k];
123449b5e25fSSatish Balay       for (j = jmin; j < jmax; j++) {
1235997a0949SHong Zhang         col       = rip[aj[j]];
1236997a0949SHong Zhang         rtmp[col] = aa[j];
12377701de02SHong Zhang         *bval++   = 0.0; /* for in-place factorization */
123849b5e25fSSatish Balay       }
12393cea4cbeSHong Zhang 
12403cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
12413cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
124249b5e25fSSatish Balay 
124391602c66SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
124449b5e25fSSatish Balay       dk = rtmp[k];
124549b5e25fSSatish Balay       i  = jl[k]; /* first row to be added to k_th row  */
124649b5e25fSSatish Balay 
1247057f5ba7SHong Zhang       while (i < k) {
124849b5e25fSSatish Balay         nexti = jl[i]; /* next row to be added to k_th row */
124949b5e25fSSatish Balay 
1250fff829cfSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
125149b5e25fSSatish Balay         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1252997a0949SHong Zhang         uikdi   = -ba[ili]*ba[bi[i]]; /* diagonal(k) */
125349b5e25fSSatish Balay         dk     += uikdi*ba[ili];
1254658e7b3eSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
125549b5e25fSSatish Balay 
1256997a0949SHong Zhang         /* add multiple of row i to k-th row */
125749b5e25fSSatish Balay         jmin = ili + 1; jmax = bi[i+1];
125849b5e25fSSatish Balay         if (jmin < jmax) {
125949b5e25fSSatish Balay           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1260dc0b31edSSatish Balay           ierr = PetscLogFlops(2.0*(jmax-jmin));CHKERRQ(ierr);
1261187a9f4bSHong Zhang 
1262fff829cfSHong Zhang           /* update il and jl for row i */
1263fff829cfSHong Zhang           il[i] = jmin;
1264fff829cfSHong Zhang           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
126549b5e25fSSatish Balay         }
1266ab27746eSHong Zhang         i = nexti;
126749b5e25fSSatish Balay       }
126849b5e25fSSatish Balay 
12693cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
12703cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
12713cea4cbeSHong Zhang       rs   = 0.0;
1272997a0949SHong Zhang       jmin = bi[k]+1;
1273997a0949SHong Zhang       nz   = bi[k+1] - jmin;
1274997a0949SHong Zhang       if (nz) {
1275997a0949SHong Zhang         bcol = bj + jmin;
1276997a0949SHong Zhang         while (nz--) {
127789f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
127889f845c8SHong Zhang           bcol++;
1279997a0949SHong Zhang         }
1280997a0949SHong Zhang       }
12813cea4cbeSHong Zhang 
12823cea4cbeSHong Zhang       sctx.rs = rs;
12833cea4cbeSHong Zhang       sctx.pv = dk;
12844cccfbddSHong Zhang       ierr    = MatPivotCheck(C,A,info,&sctx,k);CHKERRQ(ierr);
128507b50cabSHong Zhang       if (sctx.newshift) break;    /* sctx.shift_amount is updated */
12860e95ead3SHong Zhang       dk = sctx.pv;
128749b5e25fSSatish Balay 
1288997a0949SHong Zhang       /* copy data into U(k,:) */
128998a8e78dSHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
1290fff829cfSHong Zhang       jmin      = bi[k]+1; jmax = bi[k+1];
129149b5e25fSSatish Balay       if (jmin < jmax) {
129249b5e25fSSatish Balay         for (j=jmin; j<jmax; j++) {
1293997a0949SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
129449b5e25fSSatish Balay         }
129598a8e78dSHong Zhang         /* add the k-th row into il and jl */
129649b5e25fSSatish Balay         il[k] = jmin;
129798a8e78dSHong Zhang         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
129849b5e25fSSatish Balay       }
129949b5e25fSSatish Balay     }
130007b50cabSHong Zhang   } while (sctx.newshift);
1301d8c74875SBarry Smith   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
130298a8e78dSHong Zhang   if (a->permute) {ierr = PetscFree(aa);CHKERRQ(ierr);}
130349b5e25fSSatish Balay 
130449b5e25fSSatish Balay   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
130526fbe8dcSKarl Rupp 
13060a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_inplace;
13074f79d315SHong Zhang   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
13080a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace;
13090a3c351aSHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_inplace;
13100a3c351aSHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_inplace;
131149b5e25fSSatish Balay   C->assembled           = PETSC_TRUE;
13125c0bcdfcSHong Zhang   C->preallocated        = PETSC_TRUE;
131326fbe8dcSKarl Rupp 
1314d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
13153cea4cbeSHong Zhang   if (sctx.nshift) {
1316f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
131757622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1318f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
131957622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1320997a0949SHong Zhang     }
1321997a0949SHong Zhang   }
132249b5e25fSSatish Balay   PetscFunctionReturn(0);
132349b5e25fSSatish Balay }
132449b5e25fSSatish Balay 
1325671cb588SHong Zhang /*
132680d3d5a7SHong Zhang   Version for when blocks are 1 by 1 Using natural ordering under new datastructure
132780d3d5a7SHong Zhang   Modified from MatCholeskyFactorNumeric_SeqAIJ()
1328671cb588SHong Zhang */
13294a2ae208SSatish Balay #undef __FUNCT__
13304a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering"
1331d595f711SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
1332d595f711SHong Zhang {
1333d595f711SHong Zhang   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data;
13347b056e98SHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)B->data;
1335d595f711SHong Zhang   PetscErrorCode ierr;
1336d595f711SHong Zhang   PetscInt       i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bdiag=b->diag,*bjtmp;
1337545dd064SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*ajtmp;
1338d595f711SHong Zhang   PetscInt       k,jmin,jmax,*c2r,*il,col,nexti,ili,nz;
1339d595f711SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1340d90ac83dSHong Zhang   FactorShiftCtx sctx;
1341d595f711SHong Zhang   PetscReal      rs;
1342d595f711SHong Zhang   MatScalar      d,*v;
1343d595f711SHong Zhang 
1344d595f711SHong Zhang   PetscFunctionBegin;
1345dcca6d9dSJed Brown   ierr = PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&c2r);CHKERRQ(ierr);
1346545dd064SHong Zhang 
1347d595f711SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
1348d90ac83dSHong Zhang   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
1349d595f711SHong Zhang 
1350f4db908eSBarry Smith   if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
1351d595f711SHong Zhang     sctx.shift_top = info->zeropivot;
135226fbe8dcSKarl Rupp 
1353545dd064SHong Zhang     ierr = PetscMemzero(rtmp,mbs*sizeof(MatScalar));CHKERRQ(ierr);
135426fbe8dcSKarl Rupp 
1355d595f711SHong Zhang     for (i=0; i<mbs; i++) {
1356d595f711SHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
1357d595f711SHong Zhang       d        = (aa)[a->diag[i]];
1358545dd064SHong Zhang       rtmp[i] += -PetscRealPart(d);  /* diagonal entry */
1359545dd064SHong Zhang       ajtmp    = aj + ai[i] + 1;     /* exclude diagonal */
1360545dd064SHong Zhang       v        = aa + ai[i] + 1;
1361545dd064SHong Zhang       nz       = ai[i+1] - ai[i] - 1;
1362545dd064SHong Zhang       for (j=0; j<nz; j++) {
1363545dd064SHong Zhang         rtmp[i]        += PetscAbsScalar(v[j]);
1364545dd064SHong Zhang         rtmp[ajtmp[j]] += PetscAbsScalar(v[j]);
1365545dd064SHong Zhang       }
13660838c725SBarry Smith       if (PetscRealPart(rtmp[i]) > sctx.shift_top) sctx.shift_top = PetscRealPart(rtmp[i]);
1367d595f711SHong Zhang     }
1368d595f711SHong Zhang     sctx.shift_top *= 1.1;
1369d595f711SHong Zhang     sctx.nshift_max = 5;
1370d595f711SHong Zhang     sctx.shift_lo   = 0.;
1371d595f711SHong Zhang     sctx.shift_hi   = 1.;
1372d595f711SHong Zhang   }
1373d595f711SHong Zhang 
1374d595f711SHong Zhang   /* allocate working arrays
1375d595f711SHong Zhang      c2r: linked list, keep track of pivot rows for a given column. c2r[col]: head of the list for a given col
1376d595f711SHong 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
1377d595f711SHong Zhang   */
1378d595f711SHong Zhang   do {
137907b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
1380d595f711SHong Zhang 
1381d595f711SHong Zhang     for (i=0; i<mbs; i++) c2r[i] = mbs;
13829e95ef84SSatish Balay     if (mbs) il[0] = 0;
1383d595f711SHong Zhang 
1384d595f711SHong Zhang     for (k = 0; k<mbs; k++) {
1385d595f711SHong Zhang       /* zero rtmp */
1386d595f711SHong Zhang       nz    = bi[k+1] - bi[k];
1387d595f711SHong Zhang       bjtmp = bj + bi[k];
13887b056e98SHong Zhang       for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
1389d595f711SHong Zhang 
1390d595f711SHong Zhang       /* load in initial unfactored row */
1391d595f711SHong Zhang       bval = ba + bi[k];
1392d595f711SHong Zhang       jmin = ai[k]; jmax = ai[k+1];
1393d595f711SHong Zhang       for (j = jmin; j < jmax; j++) {
1394d595f711SHong Zhang         col       = aj[j];
1395d595f711SHong Zhang         rtmp[col] = aa[j];
1396d595f711SHong Zhang         *bval++   = 0.0; /* for in-place factorization */
1397d595f711SHong Zhang       }
1398d595f711SHong Zhang       /* shift the diagonal of the matrix: ZeropivotApply() */
1399d595f711SHong Zhang       rtmp[k] += sctx.shift_amount;  /* shift the diagonal of the matrix */
1400d595f711SHong Zhang 
1401d595f711SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1402d595f711SHong Zhang       dk = rtmp[k];
1403d595f711SHong Zhang       i  = c2r[k]; /* first row to be added to k_th row  */
1404d595f711SHong Zhang 
1405d595f711SHong Zhang       while (i < k) {
1406d595f711SHong Zhang         nexti = c2r[i]; /* next row to be added to k_th row */
1407d595f711SHong Zhang 
1408d595f711SHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
1409d595f711SHong Zhang         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1410d595f711SHong Zhang         uikdi   = -ba[ili]*ba[bdiag[i]]; /* diagonal(k) */
1411d595f711SHong Zhang         dk     += uikdi*ba[ili]; /* update diag[k] */
1412d595f711SHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
1413d595f711SHong Zhang 
1414d595f711SHong Zhang         /* add multiple of row i to k-th row */
1415d595f711SHong Zhang         jmin = ili + 1; jmax = bi[i+1];
1416d595f711SHong Zhang         if (jmin < jmax) {
1417d595f711SHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1418d595f711SHong Zhang           /* update il and c2r for row i */
1419d595f711SHong Zhang           il[i] = jmin;
1420d595f711SHong Zhang           j     = bj[jmin]; c2r[i] = c2r[j]; c2r[j] = i;
1421d595f711SHong Zhang         }
1422d595f711SHong Zhang         i = nexti;
1423d595f711SHong Zhang       }
1424d595f711SHong Zhang 
1425d595f711SHong Zhang       /* copy data into U(k,:) */
1426d595f711SHong Zhang       rs   = 0.0;
1427d595f711SHong Zhang       jmin = bi[k]; jmax = bi[k+1]-1;
1428d595f711SHong Zhang       if (jmin < jmax) {
1429d595f711SHong Zhang         for (j=jmin; j<jmax; j++) {
1430d595f711SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rs += PetscAbsScalar(ba[j]);
1431d595f711SHong Zhang         }
1432d595f711SHong Zhang         /* add the k-th row into il and c2r */
1433d595f711SHong Zhang         il[k] = jmin;
1434d595f711SHong Zhang         i     = bj[jmin]; c2r[k] = c2r[i]; c2r[i] = k;
1435d595f711SHong Zhang       }
1436d595f711SHong Zhang 
1437d595f711SHong Zhang       sctx.rs = rs;
1438d595f711SHong Zhang       sctx.pv = dk;
14394cccfbddSHong Zhang       ierr    = MatPivotCheck(B,A,info,&sctx,k);CHKERRQ(ierr);
144007b50cabSHong Zhang       if (sctx.newshift) break;
1441d595f711SHong Zhang       dk = sctx.pv;
1442d595f711SHong Zhang 
1443d595f711SHong Zhang       ba[bdiag[k]] = 1.0/dk; /* U(k,k) */
1444d595f711SHong Zhang     }
144507b50cabSHong Zhang   } while (sctx.newshift);
1446d595f711SHong Zhang 
1447d595f711SHong Zhang   ierr = PetscFree3(rtmp,il,c2r);CHKERRQ(ierr);
1448d595f711SHong Zhang 
1449d595f711SHong Zhang   B->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering;
14509dcdb97aSHong Zhang   B->ops->solves         = MatSolves_SeqSBAIJ_1;
1451d595f711SHong Zhang   B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering;
145280d3d5a7SHong Zhang   B->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
145380d3d5a7SHong Zhang   B->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;
1454d595f711SHong Zhang 
14557b056e98SHong Zhang   B->assembled    = PETSC_TRUE;
14567b056e98SHong Zhang   B->preallocated = PETSC_TRUE;
145726fbe8dcSKarl Rupp 
14587b056e98SHong Zhang   ierr = PetscLogFlops(B->rmap->n);CHKERRQ(ierr);
1459d595f711SHong Zhang 
1460d595f711SHong Zhang   /* MatPivotView() */
1461d595f711SHong Zhang   if (sctx.nshift) {
1462f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
146357622a8eSBarry Smith       ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %g, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,(double)sctx.shift_amount,(double)sctx.shift_fraction,(double)sctx.shift_top);CHKERRQ(ierr);
1464f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
146557622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1466f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) {
146757622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shift_inblocks applied %D, each shift_amount %g\n",sctx.nshift,(double)info->shiftamount);CHKERRQ(ierr);
1468d595f711SHong Zhang     }
1469d595f711SHong Zhang   }
1470d595f711SHong Zhang   PetscFunctionReturn(0);
1471d595f711SHong Zhang }
1472d595f711SHong Zhang 
1473d595f711SHong Zhang #undef __FUNCT__
1474d595f711SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace"
1475d595f711SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
1476671cb588SHong Zhang {
1477671cb588SHong Zhang   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ*)C->data;
1478dfbe8321SBarry Smith   PetscErrorCode ierr;
147913f74950SBarry Smith   PetscInt       i,j,mbs = a->mbs;
148013f74950SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
14813cea4cbeSHong Zhang   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
1482653a6975SHong Zhang   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
148314d2772eSHong Zhang   PetscReal      rs;
14840e95ead3SHong Zhang   FactorShiftCtx sctx;
1485653a6975SHong Zhang 
1486653a6975SHong Zhang   PetscFunctionBegin;
14870e95ead3SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
14880e95ead3SHong Zhang   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
14890e95ead3SHong Zhang 
1490653a6975SHong Zhang   /* initialization */
1491653a6975SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
1492653a6975SHong Zhang      window U(0:k, k:mbs-1).
1493653a6975SHong Zhang      jl:    list of rows to be added to uneliminated rows
1494653a6975SHong Zhang             i>= k: jl(i) is the first row to be added to row i
1495653a6975SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
1496653a6975SHong Zhang             jl(i) = mbs indicates the end of a list
1497653a6975SHong Zhang      il(i): points to the first nonzero element in U(i,k:mbs-1)
1498653a6975SHong Zhang   */
1499785e854fSJed Brown   ierr = PetscMalloc1(mbs,&rtmp);CHKERRQ(ierr);
1500dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
1501b00f7748SHong Zhang 
1502b00f7748SHong Zhang   do {
150307b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
1504*6df5ee2eSHong Zhang     il[0] = 0;
1505653a6975SHong Zhang     for (i=0; i<mbs; i++) {
1506*6df5ee2eSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs;
1507653a6975SHong Zhang     }
1508653a6975SHong Zhang 
1509997a0949SHong Zhang     for (k = 0; k<mbs; k++) {
1510653a6975SHong Zhang       /*initialize k-th row with elements nonzero in row perm(k) of A */
1511653a6975SHong Zhang       nz   = ai[k+1] - ai[k];
1512653a6975SHong Zhang       acol = aj + ai[k];
1513653a6975SHong Zhang       aval = aa + ai[k];
1514653a6975SHong Zhang       bval = ba + bi[k];
1515653a6975SHong Zhang       while (nz--) {
1516653a6975SHong Zhang         rtmp[*acol++] = *aval++;
1517653a6975SHong Zhang         *bval++       = 0.0; /* for in-place factorization */
1518653a6975SHong Zhang       }
15193cea4cbeSHong Zhang 
15203cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
15213cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1522653a6975SHong Zhang 
1523653a6975SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1524653a6975SHong Zhang       dk = rtmp[k];
1525653a6975SHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
1526653a6975SHong Zhang 
1527653a6975SHong Zhang       while (i < k) {
1528653a6975SHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
1529653a6975SHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
1530653a6975SHong Zhang         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1531653a6975SHong Zhang         uikdi   = -ba[ili]*ba[bi[i]];
1532653a6975SHong Zhang         dk     += uikdi*ba[ili];
1533653a6975SHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
1534653a6975SHong Zhang 
1535653a6975SHong Zhang         /* add multiple of row i to k-th row ... */
1536653a6975SHong Zhang         jmin = ili + 1;
1537653a6975SHong Zhang         nz   = bi[i+1] - jmin;
1538653a6975SHong Zhang         if (nz > 0) {
1539653a6975SHong Zhang           bcol = bj + jmin;
1540653a6975SHong Zhang           bval = ba + jmin;
1541dc0b31edSSatish Balay           ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
1542653a6975SHong Zhang           while (nz--) rtmp[*bcol++] += uikdi*(*bval++);
1543187a9f4bSHong Zhang 
1544997a0949SHong Zhang           /* update il and jl for i-th row */
1545997a0949SHong Zhang           il[i] = jmin;
1546997a0949SHong Zhang           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1547653a6975SHong Zhang         }
1548653a6975SHong Zhang         i = nexti;
1549653a6975SHong Zhang       }
1550653a6975SHong Zhang 
15513cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
15523cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
15533cea4cbeSHong Zhang       rs   = 0.0;
155421c26570Svictorle       jmin = bi[k]+1;
155521c26570Svictorle       nz   = bi[k+1] - jmin;
155621c26570Svictorle       if (nz) {
155721c26570Svictorle         bcol = bj + jmin;
155821c26570Svictorle         while (nz--) {
155989f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
156089f845c8SHong Zhang           bcol++;
156121c26570Svictorle         }
156221c26570Svictorle       }
15633cea4cbeSHong Zhang 
15643cea4cbeSHong Zhang       sctx.rs = rs;
15653cea4cbeSHong Zhang       sctx.pv = dk;
15664cccfbddSHong Zhang       ierr    = MatPivotCheck(C,A,info,&sctx,k);CHKERRQ(ierr);
156707b50cabSHong Zhang       if (sctx.newshift) break;    /* sctx.shift_amount is updated */
15680e95ead3SHong Zhang       dk = sctx.pv;
1569653a6975SHong Zhang 
1570997a0949SHong Zhang       /* copy data into U(k,:) */
1571653a6975SHong Zhang       ba[bi[k]] = 1.0/dk;
1572653a6975SHong Zhang       jmin      = bi[k]+1;
1573653a6975SHong Zhang       nz        = bi[k+1] - jmin;
1574653a6975SHong Zhang       if (nz) {
1575653a6975SHong Zhang         bcol = bj + jmin;
1576653a6975SHong Zhang         bval = ba + jmin;
1577653a6975SHong Zhang         while (nz--) {
1578653a6975SHong Zhang           *bval++       = rtmp[*bcol];
1579653a6975SHong Zhang           rtmp[*bcol++] = 0.0;
1580653a6975SHong Zhang         }
1581997a0949SHong Zhang         /* add k-th row into il and jl */
1582653a6975SHong Zhang         il[k] = jmin;
1583997a0949SHong Zhang         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1584653a6975SHong Zhang       }
1585b00f7748SHong Zhang     } /* end of for (k = 0; k<mbs; k++) */
158607b50cabSHong Zhang   } while (sctx.newshift);
1587653a6975SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1588d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
1589653a6975SHong Zhang 
15900a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15914f79d315SHong Zhang   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
15920a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15930a3c351aSHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15940a3c351aSHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1595db4efbfdSBarry Smith 
1596653a6975SHong Zhang   C->assembled    = PETSC_TRUE;
1597653a6975SHong Zhang   C->preallocated = PETSC_TRUE;
159826fbe8dcSKarl Rupp 
1599d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
16003cea4cbeSHong Zhang   if (sctx.nshift) {
1601f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
160257622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1603f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
160457622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1605b00f7748SHong Zhang     }
1606fb10cecfSBarry Smith   }
1607653a6975SHong Zhang   PetscFunctionReturn(0);
1608653a6975SHong Zhang }
1609653a6975SHong Zhang 
1610653a6975SHong Zhang #undef __FUNCT__
16114a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ"
16120481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,const MatFactorInfo *info)
161349b5e25fSSatish Balay {
1614dfbe8321SBarry Smith   PetscErrorCode ierr;
161549b5e25fSSatish Balay   Mat            C;
161649b5e25fSSatish Balay 
161749b5e25fSSatish Balay   PetscFunctionBegin;
1618719d5645SBarry Smith   ierr = MatGetFactor(A,"petsc",MAT_FACTOR_CHOLESKY,&C);CHKERRQ(ierr);
1619719d5645SBarry Smith   ierr = MatCholeskyFactorSymbolic(C,A,perm,info);CHKERRQ(ierr);
1620719d5645SBarry Smith   ierr = MatCholeskyFactorNumeric(C,A,info);CHKERRQ(ierr);
162126fbe8dcSKarl Rupp 
1622db4efbfdSBarry Smith   A->ops->solve          = C->ops->solve;
1623db4efbfdSBarry Smith   A->ops->solvetranspose = C->ops->solvetranspose;
162426fbe8dcSKarl Rupp 
162528be2f97SBarry Smith   ierr = MatHeaderMerge(A,&C);CHKERRQ(ierr);
162649b5e25fSSatish Balay   PetscFunctionReturn(0);
162749b5e25fSSatish Balay }
162849b5e25fSSatish Balay 
162949b5e25fSSatish Balay 
1630