xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact.c (revision 84648c2d64fd70938d2f71c337cda3739c77ab03)
1c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
2c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h>
3af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h>
4c6db04a5SJed Brown #include <petscis.h>
549b5e25fSSatish Balay 
6d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetInertia_SeqSBAIJ(Mat F, PetscInt *nneg, PetscInt *nzero, PetscInt *npos)
7d71ae5a4SJacob Faibussowitsch {
84ff4e45dSHong Zhang   Mat_SeqSBAIJ *fact = (Mat_SeqSBAIJ *)F->data;
94ff4e45dSHong Zhang   MatScalar    *dd   = fact->a;
104ff4e45dSHong Zhang   PetscInt      mbs = fact->mbs, bs = F->rmap->bs, i, nneg_tmp, npos_tmp, *fi = fact->diag;
115f9f512dSHong Zhang 
125f9f512dSHong Zhang   PetscFunctionBegin;
135f80ce2aSJacob Faibussowitsch   PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for bs: %" PetscInt_FMT " >1 yet", bs);
144ff4e45dSHong Zhang 
159371c9d4SSatish Balay   nneg_tmp = 0;
169371c9d4SSatish Balay   npos_tmp = 0;
174b8b6f3bSPierre Jolivet   if (fi) {
18eeeff2ecSHong Zhang     for (i = 0; i < mbs; i++) {
1926fbe8dcSKarl Rupp       if (PetscRealPart(dd[*fi]) > 0.0) npos_tmp++;
204ff4e45dSHong Zhang       else if (PetscRealPart(dd[*fi]) < 0.0) nneg_tmp++;
21eeeff2ecSHong Zhang       fi++;
223e0d88b5SBarry Smith     }
234b8b6f3bSPierre Jolivet   } else {
244b8b6f3bSPierre Jolivet     for (i = 0; i < mbs; i++) {
254b8b6f3bSPierre Jolivet       if (PetscRealPart(dd[fact->i[i]]) > 0.0) npos_tmp++;
264b8b6f3bSPierre Jolivet       else if (PetscRealPart(dd[fact->i[i]]) < 0.0) nneg_tmp++;
274b8b6f3bSPierre Jolivet     }
284b8b6f3bSPierre Jolivet   }
294ff4e45dSHong Zhang   if (nneg) *nneg = nneg_tmp;
30eeeff2ecSHong Zhang   if (npos) *npos = npos_tmp;
314ff4e45dSHong Zhang   if (nzero) *nzero = mbs - nneg_tmp - npos_tmp;
323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
335f9f512dSHong Zhang }
345f9f512dSHong Zhang 
355f9f512dSHong Zhang /*
365f9f512dSHong Zhang   Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP.
3710c27e3dSHong Zhang   Use Modified Sparse Row (MSR) storage for u and ju. See page 85, "Iterative Methods ..." by Saad.
385f9f512dSHong Zhang */
3966976f2fSJacob Faibussowitsch static PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(Mat F, Mat A, IS perm, const MatFactorInfo *info)
40d71ae5a4SJacob Faibussowitsch {
4110c27e3dSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ *)A->data, *b;
425d0c19d7SBarry Smith   const PetscInt *rip, *ai, *aj;
435d0c19d7SBarry Smith   PetscInt        i, mbs = a->mbs, *jutmp, bs = A->rmap->bs, bs2 = a->bs2;
4410c27e3dSHong Zhang   PetscInt        m, reallocs = 0, prow;
4510c27e3dSHong Zhang   PetscInt       *jl, *q, jmin, jmax, juidx, nzk, qm, *iu, *ju, k, j, vj, umax, maxadd;
4610c27e3dSHong Zhang   PetscReal       f = info->fill;
47ace3abfcSBarry Smith   PetscBool       perm_identity;
4810c27e3dSHong Zhang 
4910c27e3dSHong Zhang   PetscFunctionBegin;
5010c27e3dSHong Zhang   /* check whether perm is the identity mapping */
519566063dSJacob Faibussowitsch   PetscCall(ISIdentity(perm, &perm_identity));
529566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(perm, &rip));
5310c27e3dSHong Zhang 
5410c27e3dSHong Zhang   if (perm_identity) { /* without permutation */
5510c27e3dSHong Zhang     a->permute = PETSC_FALSE;
5626fbe8dcSKarl Rupp 
579371c9d4SSatish Balay     ai = a->i;
589371c9d4SSatish Balay     aj = a->j;
5910c27e3dSHong Zhang   } else { /* non-trivial permutation */
6010c27e3dSHong Zhang     a->permute = PETSC_TRUE;
6126fbe8dcSKarl Rupp 
629566063dSJacob Faibussowitsch     PetscCall(MatReorderingSeqSBAIJ(A, perm));
6326fbe8dcSKarl Rupp 
649371c9d4SSatish Balay     ai = a->inew;
659371c9d4SSatish Balay     aj = a->jnew;
6610c27e3dSHong Zhang   }
6710c27e3dSHong Zhang 
6810c27e3dSHong Zhang   /* initialization */
699f0612e4SBarry Smith   PetscCall(PetscShmgetAllocateArray(mbs + 1, sizeof(PetscInt), (void **)&iu));
709371c9d4SSatish Balay   umax = (PetscInt)(f * ai[mbs] + 1);
719371c9d4SSatish Balay   umax += mbs + 1;
729f0612e4SBarry Smith   PetscCall(PetscShmgetAllocateArray(umax, sizeof(PetscInt), (void **)&ju));
7310c27e3dSHong Zhang   iu[0] = mbs + 1;
7410c27e3dSHong Zhang   juidx = mbs + 1; /* index for ju */
75d8c74875SBarry Smith   /* jl linked list for pivot row -- linked list for col index */
769566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(mbs, &jl, mbs, &q));
7710c27e3dSHong Zhang   for (i = 0; i < mbs; i++) {
7810c27e3dSHong Zhang     jl[i] = mbs;
7910c27e3dSHong Zhang     q[i]  = 0;
8010c27e3dSHong Zhang   }
8110c27e3dSHong Zhang 
8210c27e3dSHong Zhang   /* for each row k */
8310c27e3dSHong Zhang   for (k = 0; k < mbs; k++) {
8410c27e3dSHong Zhang     for (i = 0; i < mbs; i++) q[i] = 0; /* to be removed! */
8510c27e3dSHong Zhang     nzk  = 0;                           /* num. of nz blocks in k-th block row with diagonal block excluded */
8610c27e3dSHong Zhang     q[k] = mbs;
8710c27e3dSHong Zhang     /* initialize nonzero structure of k-th row to row rip[k] of A */
8810c27e3dSHong Zhang     jmin = ai[rip[k]] + 1; /* exclude diag[k] */
8910c27e3dSHong Zhang     jmax = ai[rip[k] + 1];
9010c27e3dSHong Zhang     for (j = jmin; j < jmax; j++) {
9110c27e3dSHong Zhang       vj = rip[aj[j]]; /* col. value */
9210c27e3dSHong Zhang       if (vj > k) {
9310c27e3dSHong Zhang         qm = k;
9410c27e3dSHong Zhang         do {
959371c9d4SSatish Balay           m  = qm;
969371c9d4SSatish Balay           qm = q[m];
9710c27e3dSHong Zhang         } while (qm < vj);
985f80ce2aSJacob Faibussowitsch         PetscCheck(qm != vj, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Duplicate entry in A");
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 */
1129371c9d4SSatish Balay       jmin = iu[prow] + 1;
1139371c9d4SSatish Balay       jmax = iu[prow + 1];
11410c27e3dSHong Zhang       qm   = k;
11510c27e3dSHong Zhang       for (j = jmin; j < jmax; j++) {
11610c27e3dSHong Zhang         vj = ju[j];
11710c27e3dSHong Zhang         do {
1189371c9d4SSatish Balay           m  = qm;
1199371c9d4SSatish Balay           qm = q[m];
12010c27e3dSHong Zhang         } while (qm < vj);
12110c27e3dSHong Zhang         if (qm != vj) {
1229371c9d4SSatish Balay           nzk++;
1239371c9d4SSatish Balay           q[m]  = vj;
1249371c9d4SSatish Balay           q[vj] = qm;
1259371c9d4SSatish Balay           qm    = vj;
12610c27e3dSHong Zhang         }
12710c27e3dSHong Zhang       }
12810c27e3dSHong Zhang       prow = jl[prow]; /* next pivot row */
12910c27e3dSHong Zhang     }
13010c27e3dSHong Zhang 
13110c27e3dSHong Zhang     /* add k to row list for first nonzero element in k-th row */
13210c27e3dSHong Zhang     if (nzk > 0) {
13310c27e3dSHong Zhang       i     = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */
1349371c9d4SSatish Balay       jl[k] = jl[i];
1359371c9d4SSatish Balay       jl[i] = k;
13610c27e3dSHong Zhang     }
13710c27e3dSHong Zhang     iu[k + 1] = iu[k] + nzk;
13810c27e3dSHong Zhang 
13910c27e3dSHong Zhang     /* allocate more space to ju if needed */
14010c27e3dSHong Zhang     if (iu[k + 1] > umax) {
14110c27e3dSHong Zhang       /* estimate how much additional space we will need */
14210c27e3dSHong Zhang       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
14310c27e3dSHong Zhang       /* just double the memory each time */
14410c27e3dSHong Zhang       maxadd = umax;
14510c27e3dSHong Zhang       if (maxadd < nzk) maxadd = (mbs - k) * (nzk + 1) / 2;
14610c27e3dSHong Zhang       umax += maxadd;
14710c27e3dSHong Zhang 
14810c27e3dSHong Zhang       /* allocate a longer ju */
1499f0612e4SBarry Smith       PetscCall(PetscShmgetAllocateArray(umax, sizeof(PetscInt), (void **)&jutmp));
1509566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(jutmp, ju, iu[k]));
1519f0612e4SBarry Smith       PetscCall(PetscShmgetDeallocateArray((void **)&ju));
15210c27e3dSHong Zhang       ju = jutmp;
15310c27e3dSHong Zhang       reallocs++; /* count how many times we realloc */
15410c27e3dSHong Zhang     }
15510c27e3dSHong Zhang 
15610c27e3dSHong Zhang     /* save nonzero structure of k-th row in ju */
15710c27e3dSHong Zhang     i = k;
15810c27e3dSHong Zhang     while (nzk--) {
15910c27e3dSHong Zhang       i           = q[i];
16010c27e3dSHong Zhang       ju[juidx++] = i;
16110c27e3dSHong Zhang     }
16210c27e3dSHong Zhang   }
16310c27e3dSHong Zhang 
1646cf91177SBarry Smith #if defined(PETSC_USE_INFO)
16510c27e3dSHong Zhang   if (ai[mbs] != 0) {
16610c27e3dSHong Zhang     PetscReal af = ((PetscReal)iu[mbs]) / ((PetscReal)ai[mbs]);
1679566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n", reallocs, (double)f, (double)af));
1689566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Run with -pc_factor_fill %g or use \n", (double)af));
1699566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "PCFactorSetFill(pc,%g);\n", (double)af));
1709566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "for best performance.\n"));
17110c27e3dSHong Zhang   } else {
1729566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Empty matrix\n"));
17310c27e3dSHong Zhang   }
17463ba0a88SBarry Smith #endif
17510c27e3dSHong Zhang 
1769566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(perm, &rip));
1779566063dSJacob Faibussowitsch   PetscCall(PetscFree2(jl, q));
17810c27e3dSHong Zhang 
17910c27e3dSHong Zhang   /* put together the new matrix */
1809566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(F, bs, MAT_SKIP_ALLOCATION, NULL));
18110c27e3dSHong Zhang 
18257508eceSPierre Jolivet   b          = (Mat_SeqSBAIJ *)F->data;
183e6b907acSBarry Smith   b->free_ij = PETSC_TRUE;
1849f0612e4SBarry Smith   PetscCall(PetscShmgetAllocateArray((iu[mbs] + 1) * bs2, sizeof(PetscScalar), (void **)&b->a));
1859f0612e4SBarry Smith   b->free_a = PETSC_TRUE;
18610c27e3dSHong Zhang   b->j      = ju;
18710c27e3dSHong Zhang   b->i      = iu;
188f4259b30SLisandro Dalcin   b->diag   = NULL;
189f4259b30SLisandro Dalcin   b->ilen   = NULL;
190f4259b30SLisandro Dalcin   b->imax   = NULL;
19110c27e3dSHong Zhang   b->row    = perm;
19226fbe8dcSKarl Rupp 
19310c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
19426fbe8dcSKarl Rupp 
1959566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)perm));
19626fbe8dcSKarl Rupp 
19710c27e3dSHong Zhang   b->icol = perm;
1989566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)perm));
1999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs * mbs + bs, &b->solve_work));
20010c27e3dSHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.
20110c27e3dSHong Zhang      Allocate idnew, solve_work, new a, new j */
20210c27e3dSHong Zhang   b->maxnz = b->nz = iu[mbs];
20310c27e3dSHong Zhang 
20457508eceSPierre Jolivet   F->info.factor_mallocs   = reallocs;
20557508eceSPierre Jolivet   F->info.fill_ratio_given = f;
20610c27e3dSHong Zhang   if (ai[mbs] != 0) {
20757508eceSPierre Jolivet     F->info.fill_ratio_needed = ((PetscReal)iu[mbs]) / ((PetscReal)ai[mbs]);
20810c27e3dSHong Zhang   } else {
20957508eceSPierre Jolivet     F->info.fill_ratio_needed = 0.0;
21010c27e3dSHong Zhang   }
2119566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetNumericFactorization_inplace(F, perm_identity));
2123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21310c27e3dSHong Zhang }
21410c27e3dSHong Zhang /*
21510c27e3dSHong Zhang     Symbolic U^T*D*U factorization for SBAIJ format.
21680d3d5a7SHong Zhang     See MatICCFactorSymbolic_SeqAIJ() for description of its data structure.
21710c27e3dSHong Zhang */
218c6db04a5SJed Brown #include <petscbt.h>
219c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
220d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat fact, Mat A, IS perm, const MatFactorInfo *info)
221d71ae5a4SJacob Faibussowitsch {
22298a8e78dSHong Zhang   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ *)A->data;
22398a8e78dSHong Zhang   Mat_SeqSBAIJ      *b;
224ace3abfcSBarry Smith   PetscBool          perm_identity, missing;
22598a8e78dSHong Zhang   PetscReal          fill = info->fill;
2267b056e98SHong Zhang   const PetscInt    *rip, *ai = a->i, *aj = a->j;
2279186b105SHong Zhang   PetscInt           i, mbs = a->mbs, bs = A->rmap->bs, reallocs = 0, prow;
22898a8e78dSHong Zhang   PetscInt          *jl, jmin, jmax, nzk, *ui, k, j, *il, nextprow;
229c6d0d4f0SHong Zhang   PetscInt           nlnk, *lnk, ncols, *cols, *uj, **ui_ptr, *uj_ptr, *udiag;
2300298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL, current_space = NULL;
23198a8e78dSHong Zhang   PetscBT            lnkbt;
232d595f711SHong Zhang 
233d595f711SHong Zhang   PetscFunctionBegin;
2345f80ce2aSJacob Faibussowitsch   PetscCheck(A->rmap->n == A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Must be square matrix, rows %" PetscInt_FMT " columns %" PetscInt_FMT, A->rmap->n, A->cmap->n);
2359566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal(A, &missing, &i));
2365f80ce2aSJacob Faibussowitsch   PetscCheck(!missing, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal entry %" PetscInt_FMT, i);
2376d07c2b1SHong Zhang   if (bs > 1) {
2389566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(fact, A, perm, info));
2393ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2406d07c2b1SHong Zhang   }
241d595f711SHong Zhang 
242d595f711SHong Zhang   /* check whether perm is the identity mapping */
2439566063dSJacob Faibussowitsch   PetscCall(ISIdentity(perm, &perm_identity));
2445f80ce2aSJacob Faibussowitsch   PetscCheck(perm_identity, PETSC_COMM_SELF, PETSC_ERR_SUP, "Matrix reordering is not supported for sbaij matrix. Use aij format");
245d595f711SHong Zhang   a->permute = PETSC_FALSE;
2469566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(perm, &rip));
247d595f711SHong Zhang 
248d595f711SHong Zhang   /* initialization */
2499f0612e4SBarry Smith   PetscCall(PetscShmgetAllocateArray(mbs + 1, sizeof(PetscInt), (void **)&ui));
2509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs + 1, &udiag));
251d595f711SHong Zhang   ui[0] = 0;
252d595f711SHong Zhang 
253d595f711SHong Zhang   /* jl: linked list for storing indices of the pivot rows
254d595f711SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
2559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(mbs, &ui_ptr, mbs, &il, mbs, &jl, mbs, &cols));
256d595f711SHong Zhang   for (i = 0; i < mbs; i++) {
2579371c9d4SSatish Balay     jl[i] = mbs;
2589371c9d4SSatish Balay     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;
2639566063dSJacob Faibussowitsch   PetscCall(PetscLLCreate(mbs, mbs, nlnk, lnk, lnkbt));
264d595f711SHong Zhang 
265d595f711SHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
2669566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill, ai[mbs] + 1), &free_space));
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];
2735f80ce2aSJacob Faibussowitsch     PetscCheck(ncols, PETSC_COMM_SELF, PETSC_ERR_MAT_CH_ZRPVT, "Empty row %" PetscInt_FMT " 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     }
2789566063dSJacob Faibussowitsch     PetscCall(PetscLLAdd(ncols, cols, mbs, &nlnk, lnk, lnkbt));
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) */
2919566063dSJacob Faibussowitsch       PetscCall(PetscLLAddSorted(ncols, uj_ptr, mbs, &nlnk, lnk, lnkbt));
292d595f711SHong Zhang       nzk += nlnk;
293d595f711SHong Zhang 
294d595f711SHong Zhang       /* update il and jl for prow */
295d595f711SHong Zhang       if (jmin < jmax) {
296d595f711SHong Zhang         il[prow] = jmin;
2979371c9d4SSatish Balay         j        = *uj_ptr;
2989371c9d4SSatish Balay         jl[prow] = jl[j];
2999371c9d4SSatish Balay         jl[j]    = prow;
300d595f711SHong Zhang       }
301d595f711SHong Zhang       prow = nextprow;
302d595f711SHong Zhang     }
303d595f711SHong Zhang 
304d595f711SHong Zhang     /* if free space is not available, make more free space */
305d595f711SHong Zhang     if (current_space->local_remaining < nzk) {
306d595f711SHong Zhang       i = mbs - k + 1;                                   /* num of unfactored rows */
307f91af8c7SBarry Smith       i = PetscIntMultTruncate(i, PetscMin(nzk, i - 1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
3089566063dSJacob Faibussowitsch       PetscCall(PetscFreeSpaceGet(i, &current_space));
309d595f711SHong Zhang       reallocs++;
310d595f711SHong Zhang     }
311d595f711SHong Zhang 
312d595f711SHong Zhang     /* copy data into free space, then initialize lnk */
3139566063dSJacob Faibussowitsch     PetscCall(PetscLLClean(mbs, mbs, nzk, lnk, current_space->array, lnkbt));
314d595f711SHong Zhang 
315d595f711SHong Zhang     /* add the k-th row into il and jl */
3167b056e98SHong Zhang     if (nzk > 1) {
317d595f711SHong Zhang       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
3189371c9d4SSatish Balay       jl[k] = jl[i];
3199371c9d4SSatish Balay       jl[i] = k;
320d595f711SHong Zhang       il[k] = ui[k] + 1;
321d595f711SHong Zhang     }
322d595f711SHong Zhang     ui_ptr[k] = current_space->array;
32326fbe8dcSKarl Rupp 
324d595f711SHong Zhang     current_space->array += nzk;
325d595f711SHong Zhang     current_space->local_used += nzk;
326d595f711SHong Zhang     current_space->local_remaining -= nzk;
327d595f711SHong Zhang 
328d595f711SHong Zhang     ui[k + 1] = ui[k] + nzk;
329d595f711SHong Zhang   }
330d595f711SHong Zhang 
3319566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(perm, &rip));
3329566063dSJacob Faibussowitsch   PetscCall(PetscFree4(ui_ptr, il, jl, cols));
333d595f711SHong Zhang 
334d595f711SHong Zhang   /* destroy list of free space and other temporary array(s) */
335*84648c2dSPierre Jolivet   PetscCall(PetscShmgetAllocateArray(ui[mbs], sizeof(PetscInt), (void **)&uj));
3369566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceContiguous_Cholesky(&free_space, uj, mbs, ui, udiag)); /* store matrix factor */
3379566063dSJacob Faibussowitsch   PetscCall(PetscLLDestroy(lnk, lnkbt));
338d595f711SHong Zhang 
339d595f711SHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
3409566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(fact, bs, MAT_SKIP_ALLOCATION, NULL));
341d595f711SHong Zhang 
3427b056e98SHong Zhang   b          = (Mat_SeqSBAIJ *)fact->data;
343d595f711SHong Zhang   b->free_ij = PETSC_TRUE;
344*84648c2dSPierre Jolivet   PetscCall(PetscShmgetAllocateArray(ui[mbs], sizeof(PetscScalar), (void **)&b->a));
3459f0612e4SBarry Smith   b->free_a    = PETSC_TRUE;
346d595f711SHong Zhang   b->j         = uj;
347d595f711SHong Zhang   b->i         = ui;
348c6d0d4f0SHong Zhang   b->diag      = udiag;
349c6d0d4f0SHong Zhang   b->free_diag = PETSC_TRUE;
350f4259b30SLisandro Dalcin   b->ilen      = NULL;
351f4259b30SLisandro Dalcin   b->imax      = NULL;
352d595f711SHong Zhang   b->row       = perm;
353d595f711SHong Zhang   b->icol      = perm;
35426fbe8dcSKarl Rupp 
3559566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)perm));
3569566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)perm));
35726fbe8dcSKarl Rupp 
3587b056e98SHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
35926fbe8dcSKarl Rupp 
3609566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs + 1, &b->solve_work));
36126fbe8dcSKarl Rupp 
362d595f711SHong Zhang   b->maxnz = b->nz = ui[mbs];
363d595f711SHong Zhang 
36495b5ac22SHong Zhang   fact->info.factor_mallocs   = reallocs;
36595b5ac22SHong Zhang   fact->info.fill_ratio_given = fill;
366d595f711SHong Zhang   if (ai[mbs] != 0) {
36795b5ac22SHong Zhang     fact->info.fill_ratio_needed = ((PetscReal)ui[mbs]) / ai[mbs];
368d595f711SHong Zhang   } else {
36995b5ac22SHong Zhang     fact->info.fill_ratio_needed = 0.0;
370d595f711SHong Zhang   }
37195b5ac22SHong Zhang #if defined(PETSC_USE_INFO)
37295b5ac22SHong Zhang   if (ai[mbs] != 0) {
37395b5ac22SHong Zhang     PetscReal af = fact->info.fill_ratio_needed;
3749566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n", reallocs, (double)fill, (double)af));
3759566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Run with -pc_factor_fill %g or use \n", (double)af));
3769566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "PCFactorSetFill(pc,%g) for best performance.\n", (double)af));
37795b5ac22SHong Zhang   } else {
3789566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Empty matrix\n"));
37995b5ac22SHong Zhang   }
38095b5ac22SHong Zhang #endif
381c6d0d4f0SHong Zhang   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
3823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
383d595f711SHong Zhang }
384d595f711SHong Zhang 
385d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(Mat fact, Mat A, IS perm, const MatFactorInfo *info)
386d71ae5a4SJacob Faibussowitsch {
387d595f711SHong Zhang   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ *)A->data;
388d595f711SHong Zhang   Mat_SeqSBAIJ      *b;
389ace3abfcSBarry Smith   PetscBool          perm_identity, missing;
390d595f711SHong Zhang   PetscReal          fill = info->fill;
391d595f711SHong Zhang   const PetscInt    *rip, *ai, *aj;
392d595f711SHong Zhang   PetscInt           i, mbs = a->mbs, bs = A->rmap->bs, reallocs = 0, prow, d;
393d595f711SHong Zhang   PetscInt          *jl, jmin, jmax, nzk, *ui, k, j, *il, nextprow;
394d595f711SHong Zhang   PetscInt           nlnk, *lnk, ncols, *cols, *uj, **ui_ptr, *uj_ptr;
3950298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL, current_space = NULL;
396d595f711SHong Zhang   PetscBT            lnkbt;
39749b5e25fSSatish Balay 
39849b5e25fSSatish Balay   PetscFunctionBegin;
3999566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal(A, &missing, &d));
4005f80ce2aSJacob Faibussowitsch   PetscCheck(!missing, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal entry %" PetscInt_FMT, d);
40158ebbce7SBarry Smith 
40210c27e3dSHong Zhang   /*
40310c27e3dSHong Zhang    This code originally uses Modified Sparse Row (MSR) storage
404da81f932SPierre Jolivet    (see page 85, "Iterative Methods ..." by Saad) for the output matrix B - bad choice!
40510c27e3dSHong Zhang    Then it is rewritten so the factor B takes seqsbaij format. However the associated
40610c27e3dSHong Zhang    MatCholeskyFactorNumeric_() have not been modified for the cases of bs>1 or !perm_identity,
40710c27e3dSHong Zhang    thus the original code in MSR format is still used for these cases.
40810c27e3dSHong Zhang    The code below should replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR() whenever
40910c27e3dSHong Zhang    MatCholeskyFactorNumeric_() is modified for using sbaij symbolic factor.
41010c27e3dSHong Zhang   */
411fff829cfSHong Zhang   if (bs > 1) {
4129566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(fact, A, perm, info));
4133ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
414fff829cfSHong Zhang   }
41510c27e3dSHong Zhang 
41698a8e78dSHong Zhang   /* check whether perm is the identity mapping */
4179566063dSJacob Faibussowitsch   PetscCall(ISIdentity(perm, &perm_identity));
4185f80ce2aSJacob Faibussowitsch   PetscCheck(perm_identity, PETSC_COMM_SELF, PETSC_ERR_SUP, "Matrix reordering is not supported for sbaij matrix. Use aij format");
419fff829cfSHong Zhang   a->permute = PETSC_FALSE;
4209371c9d4SSatish Balay   ai         = a->i;
4219371c9d4SSatish Balay   aj         = a->j;
4229566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(perm, &rip));
423fff829cfSHong Zhang 
424fff829cfSHong Zhang   /* initialization */
4259f0612e4SBarry Smith   PetscCall(PetscShmgetAllocateArray(mbs + 1, sizeof(PetscInt), (void **)&ui));
42698a8e78dSHong Zhang   ui[0] = 0;
42798a8e78dSHong Zhang 
42898a8e78dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
4297625dc9aSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
4309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(mbs, &ui_ptr, mbs, &il, mbs, &jl, mbs, &cols));
4317625dc9aSHong Zhang   for (i = 0; i < mbs; i++) {
4329371c9d4SSatish Balay     jl[i] = mbs;
4339371c9d4SSatish Balay     il[i] = 0;
434fff829cfSHong Zhang   }
435fff829cfSHong Zhang 
43698a8e78dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
4377625dc9aSHong Zhang   nlnk = mbs + 1;
4389566063dSJacob Faibussowitsch   PetscCall(PetscLLCreate(mbs, mbs, nlnk, lnk, lnkbt));
439fff829cfSHong Zhang 
4407625dc9aSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
4419566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill, ai[mbs] + 1), &free_space));
44298a8e78dSHong Zhang   current_space = free_space;
44398a8e78dSHong Zhang 
4447625dc9aSHong Zhang   for (k = 0; k < mbs; k++) { /* for each active row k */
44598a8e78dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
44698a8e78dSHong Zhang     nzk   = 0;
44798a8e78dSHong Zhang     ncols = ai[rip[k] + 1] - ai[rip[k]];
4487625dc9aSHong Zhang     for (j = 0; j < ncols; j++) {
4497625dc9aSHong Zhang       i       = *(aj + ai[rip[k]] + j);
4507625dc9aSHong Zhang       cols[j] = rip[i];
4517625dc9aSHong Zhang     }
4529566063dSJacob Faibussowitsch     PetscCall(PetscLLAdd(ncols, cols, mbs, &nlnk, lnk, lnkbt));
45398a8e78dSHong Zhang     nzk += nlnk;
45498a8e78dSHong Zhang 
45598a8e78dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
45698a8e78dSHong Zhang     prow = jl[k]; /* 1st pivot row */
457fff829cfSHong Zhang 
458fff829cfSHong Zhang     while (prow < k) {
459fff829cfSHong Zhang       nextprow = jl[prow];
46098a8e78dSHong Zhang       /* merge prow into k-th row */
4617625dc9aSHong Zhang       jmin   = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
46298a8e78dSHong Zhang       jmax   = ui[prow + 1];
46398a8e78dSHong Zhang       ncols  = jmax - jmin;
4647625dc9aSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
4659566063dSJacob Faibussowitsch       PetscCall(PetscLLAddSorted(ncols, uj_ptr, mbs, &nlnk, lnk, lnkbt));
46698a8e78dSHong Zhang       nzk += nlnk;
467fff829cfSHong Zhang 
46898a8e78dSHong Zhang       /* update il and jl for prow */
469fff829cfSHong Zhang       if (jmin < jmax) {
470fff829cfSHong Zhang         il[prow] = jmin;
47126fbe8dcSKarl Rupp 
4729371c9d4SSatish Balay         j        = *uj_ptr;
4739371c9d4SSatish Balay         jl[prow] = jl[j];
4749371c9d4SSatish Balay         jl[j]    = prow;
475fff829cfSHong Zhang       }
476fff829cfSHong Zhang       prow = nextprow;
477fff829cfSHong Zhang     }
478fff829cfSHong Zhang 
47998a8e78dSHong Zhang     /* if free space is not available, make more free space */
48098a8e78dSHong Zhang     if (current_space->local_remaining < nzk) {
4817625dc9aSHong Zhang       i = mbs - k + 1;                                                            /* num of unfactored rows */
482f91af8c7SBarry Smith       i = PetscMin(PetscIntMultTruncate(i, nzk), PetscIntMultTruncate(i, i - 1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
4839566063dSJacob Faibussowitsch       PetscCall(PetscFreeSpaceGet(i, &current_space));
48498a8e78dSHong Zhang       reallocs++;
48598a8e78dSHong Zhang     }
48698a8e78dSHong Zhang 
48798a8e78dSHong Zhang     /* copy data into free space, then initialize lnk */
4889566063dSJacob Faibussowitsch     PetscCall(PetscLLClean(mbs, mbs, nzk, lnk, current_space->array, lnkbt));
48998a8e78dSHong Zhang 
49098a8e78dSHong Zhang     /* add the k-th row into il and jl */
49198a8e78dSHong Zhang     if (nzk - 1 > 0) {
4927625dc9aSHong Zhang       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
4939371c9d4SSatish Balay       jl[k] = jl[i];
4949371c9d4SSatish Balay       jl[i] = k;
49598a8e78dSHong Zhang       il[k] = ui[k] + 1;
496fff829cfSHong Zhang     }
4977625dc9aSHong Zhang     ui_ptr[k] = current_space->array;
49826fbe8dcSKarl Rupp 
49998a8e78dSHong Zhang     current_space->array += nzk;
50098a8e78dSHong Zhang     current_space->local_used += nzk;
50198a8e78dSHong Zhang     current_space->local_remaining -= nzk;
502fff829cfSHong Zhang 
50398a8e78dSHong Zhang     ui[k + 1] = ui[k] + nzk;
504fff829cfSHong Zhang   }
505fff829cfSHong Zhang 
5069566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(perm, &rip));
5079566063dSJacob Faibussowitsch   PetscCall(PetscFree4(ui_ptr, il, jl, cols));
508fff829cfSHong Zhang 
50998a8e78dSHong Zhang   /* destroy list of free space and other temporary array(s) */
5109f0612e4SBarry Smith   PetscCall(PetscShmgetAllocateArray(ui[mbs] + 1, sizeof(PetscInt), (void **)&uj));
5119566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceContiguous(&free_space, uj));
5129566063dSJacob Faibussowitsch   PetscCall(PetscLLDestroy(lnk, lnkbt));
513fff829cfSHong Zhang 
51498a8e78dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
5159566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(fact, bs, MAT_SKIP_ALLOCATION, NULL));
51698a8e78dSHong Zhang 
51795b5ac22SHong Zhang   b = (Mat_SeqSBAIJ *)fact->data;
5189f0612e4SBarry Smith   PetscCall(PetscShmgetAllocateArray(ui[mbs] + 1, sizeof(PetscScalar), (void **)&b->a));
519e6b907acSBarry Smith   b->free_a  = PETSC_TRUE;
520e6b907acSBarry Smith   b->free_ij = PETSC_TRUE;
52198a8e78dSHong Zhang   b->j       = uj;
52298a8e78dSHong Zhang   b->i       = ui;
523f4259b30SLisandro Dalcin   b->diag    = NULL;
524f4259b30SLisandro Dalcin   b->ilen    = NULL;
525f4259b30SLisandro Dalcin   b->imax    = NULL;
526fff829cfSHong Zhang   b->row     = perm;
52726fbe8dcSKarl Rupp 
528fff829cfSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
52926fbe8dcSKarl Rupp 
5309566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)perm));
531fff829cfSHong Zhang   b->icol = perm;
5329566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)perm));
5339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs + 1, &b->solve_work));
5347625dc9aSHong Zhang   b->maxnz = b->nz = ui[mbs];
535fff829cfSHong Zhang 
53695b5ac22SHong Zhang   fact->info.factor_mallocs   = reallocs;
53795b5ac22SHong Zhang   fact->info.fill_ratio_given = fill;
5387625dc9aSHong Zhang   if (ai[mbs] != 0) {
53995b5ac22SHong Zhang     fact->info.fill_ratio_needed = ((PetscReal)ui[mbs]) / ai[mbs];
540fff829cfSHong Zhang   } else {
54195b5ac22SHong Zhang     fact->info.fill_ratio_needed = 0.0;
542fff829cfSHong Zhang   }
54395b5ac22SHong Zhang #if defined(PETSC_USE_INFO)
54495b5ac22SHong Zhang   if (ai[mbs] != 0) {
54595b5ac22SHong Zhang     PetscReal af = fact->info.fill_ratio_needed;
5469566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n", reallocs, (double)fill, (double)af));
5479566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Run with -pc_factor_fill %g or use \n", (double)af));
5489566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "PCFactorSetFill(pc,%g) for best performance.\n", (double)af));
54995b5ac22SHong Zhang   } else {
5509566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Empty matrix\n"));
55195b5ac22SHong Zhang   }
55295b5ac22SHong Zhang #endif
5539566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetNumericFactorization_inplace(fact, perm_identity));
5543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
555064503c5SHong Zhang }
556d595f711SHong Zhang 
557d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat C, Mat A, const MatFactorInfo *info)
558d71ae5a4SJacob Faibussowitsch {
5594c16a6a6SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ *)A->data, *b = (Mat_SeqSBAIJ *)C->data;
5604c16a6a6SHong Zhang   IS              perm = b->row;
5615d0c19d7SBarry Smith   const PetscInt *ai, *aj, *perm_ptr, mbs = a->mbs, *bi = b->i, *bj = b->j;
5625d0c19d7SBarry Smith   PetscInt        i, j;
5635d0c19d7SBarry Smith   PetscInt       *a2anew, k, k1, jmin, jmax, *jl, *il, vj, nexti, ili;
5643bc0b13bSBarry Smith   PetscInt        bs = A->rmap->bs, bs2 = a->bs2;
5654c16a6a6SHong Zhang   MatScalar      *ba = b->a, *aa, *ap, *dk, *uik;
5664c16a6a6SHong Zhang   MatScalar      *u, *diag, *rtmp, *rtmp_ptr;
56728de702eSHong Zhang   MatScalar      *work;
56813f74950SBarry Smith   PetscInt       *pivots;
5695f8bbccaSHong Zhang   PetscBool       allowzeropivot, zeropivotdetected;
5704c16a6a6SHong Zhang 
5714c16a6a6SHong Zhang   PetscFunctionBegin;
5724c16a6a6SHong Zhang   /* initialization */
5739566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(bs2 * mbs, &rtmp));
5749566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(mbs, &il, mbs, &jl));
5755f8bbccaSHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
5765f8bbccaSHong Zhang 
5776df5ee2eSHong Zhang   il[0] = 0;
5786df5ee2eSHong Zhang   for (i = 0; i < mbs; i++) jl[i] = mbs;
5796df5ee2eSHong Zhang 
5809566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(bs2, &dk, bs2, &uik, bs, &work));
5819566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs, &pivots));
5824c16a6a6SHong Zhang 
5839566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(perm, &perm_ptr));
5844c16a6a6SHong Zhang 
5854c16a6a6SHong Zhang   /* check permutation */
5864c16a6a6SHong Zhang   if (!a->permute) {
5879371c9d4SSatish Balay     ai = a->i;
5889371c9d4SSatish Balay     aj = a->j;
5899371c9d4SSatish Balay     aa = a->a;
5904c16a6a6SHong Zhang   } else {
5919371c9d4SSatish Balay     ai = a->inew;
5929371c9d4SSatish Balay     aj = a->jnew;
5939566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs2 * ai[mbs], &aa));
5949566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(aa, a->a, bs2 * ai[mbs]));
5959566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(ai[mbs], &a2anew));
5969566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(a2anew, a->a2anew, ai[mbs]));
5974c16a6a6SHong Zhang 
5984c16a6a6SHong Zhang     for (i = 0; i < mbs; i++) {
5999371c9d4SSatish Balay       jmin = ai[i];
6009371c9d4SSatish Balay       jmax = ai[i + 1];
6014c16a6a6SHong Zhang       for (j = jmin; j < jmax; j++) {
6024c16a6a6SHong Zhang         while (a2anew[j] != j) {
6039371c9d4SSatish Balay           k         = a2anew[j];
6049371c9d4SSatish Balay           a2anew[j] = a2anew[k];
6059371c9d4SSatish Balay           a2anew[k] = k;
6064c16a6a6SHong Zhang           for (k1 = 0; k1 < bs2; k1++) {
6074c16a6a6SHong Zhang             dk[k1]           = aa[k * bs2 + k1];
6084c16a6a6SHong Zhang             aa[k * bs2 + k1] = aa[j * bs2 + k1];
6094c16a6a6SHong Zhang             aa[j * bs2 + k1] = dk[k1];
6104c16a6a6SHong Zhang           }
6114c16a6a6SHong Zhang         }
6125e116b59SBarry Smith         /* transform column-oriented blocks that lie in the lower triangle to row-oriented blocks */
6134c16a6a6SHong Zhang         if (i > aj[j]) {
6144c16a6a6SHong Zhang           ap = aa + j * bs2;                       /* ptr to the beginning of j-th block of aa */
6154c16a6a6SHong Zhang           for (k = 0; k < bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
6164c16a6a6SHong Zhang           for (k = 0; k < bs; k++) {               /* j-th block of aa <- dk^T */
6174c16a6a6SHong Zhang             for (k1 = 0; k1 < bs; k1++) *ap++ = dk[k + bs * k1];
6184c16a6a6SHong Zhang           }
6194c16a6a6SHong Zhang         }
6204c16a6a6SHong Zhang       }
6214c16a6a6SHong Zhang     }
6229566063dSJacob Faibussowitsch     PetscCall(PetscFree(a2anew));
6234c16a6a6SHong Zhang   }
6244c16a6a6SHong Zhang 
6254c16a6a6SHong Zhang   /* for each row k */
6264c16a6a6SHong Zhang   for (k = 0; k < mbs; k++) {
6274c16a6a6SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
6289371c9d4SSatish Balay     jmin = ai[perm_ptr[k]];
6299371c9d4SSatish Balay     jmax = ai[perm_ptr[k] + 1];
630057f5ba7SHong Zhang 
6314c16a6a6SHong Zhang     ap = aa + jmin * bs2;
6324c16a6a6SHong Zhang     for (j = jmin; j < jmax; j++) {
6334c16a6a6SHong Zhang       vj       = perm_ptr[aj[j]]; /* block col. index */
6344c16a6a6SHong Zhang       rtmp_ptr = rtmp + vj * bs2;
6354c16a6a6SHong Zhang       for (i = 0; i < bs2; i++) *rtmp_ptr++ = *ap++;
6364c16a6a6SHong Zhang     }
6374c16a6a6SHong Zhang 
6384c16a6a6SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
6399566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dk, rtmp + k * bs2, bs2));
6404c16a6a6SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
6414c16a6a6SHong Zhang 
642057f5ba7SHong Zhang     while (i < k) {
6434c16a6a6SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
6444c16a6a6SHong Zhang 
6454c16a6a6SHong Zhang       /* compute multiplier */
6464c16a6a6SHong Zhang       ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
6474c16a6a6SHong Zhang 
6484c16a6a6SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
6494c16a6a6SHong Zhang       diag = ba + i * bs2;
6504c16a6a6SHong Zhang       u    = ba + ili * bs2;
6519566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(uik, bs2));
65296b95a6bSBarry Smith       PetscKernel_A_gets_A_minus_B_times_C(bs, uik, diag, u);
6534c16a6a6SHong Zhang 
6544c16a6a6SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
65596b95a6bSBarry Smith       PetscKernel_A_gets_A_plus_Btranspose_times_C(bs, dk, uik, u);
6569566063dSJacob Faibussowitsch       PetscCall(PetscLogFlops(4.0 * bs * bs2));
6574c16a6a6SHong Zhang 
6584c16a6a6SHong Zhang       /* update -U(i,k) */
6599566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(ba + ili * bs2, uik, bs2));
6604c16a6a6SHong Zhang 
6614c16a6a6SHong Zhang       /* add multiple of row i to k-th row ... */
6629371c9d4SSatish Balay       jmin = ili + 1;
6639371c9d4SSatish Balay       jmax = bi[i + 1];
6644c16a6a6SHong Zhang       if (jmin < jmax) {
6654c16a6a6SHong Zhang         for (j = jmin; j < jmax; j++) {
6664c16a6a6SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
6674c16a6a6SHong Zhang           rtmp_ptr = rtmp + bj[j] * bs2;
6684c16a6a6SHong Zhang           u        = ba + j * bs2;
66996b95a6bSBarry Smith           PetscKernel_A_gets_A_plus_Btranspose_times_C(bs, rtmp_ptr, uik, u);
6704c16a6a6SHong Zhang         }
6719566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(2.0 * bs * bs2 * (jmax - jmin)));
6724c16a6a6SHong Zhang 
6734c16a6a6SHong Zhang         /* ... add i to row list for next nonzero entry */
6744c16a6a6SHong Zhang         il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
6754c16a6a6SHong Zhang         j     = bj[jmin];
6769371c9d4SSatish Balay         jl[i] = jl[j];
6779371c9d4SSatish Balay         jl[j] = i; /* update jl */
6784c16a6a6SHong Zhang       }
6794c16a6a6SHong Zhang       i = nexti;
6804c16a6a6SHong Zhang     }
6814c16a6a6SHong Zhang 
6824c16a6a6SHong Zhang     /* save nonzero entries in k-th row of U ... */
6834c16a6a6SHong Zhang 
6844c16a6a6SHong Zhang     /* invert diagonal block */
6854c16a6a6SHong Zhang     diag = ba + k * bs2;
6869566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(diag, dk, bs2));
6875f8bbccaSHong Zhang 
6889566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A(bs, diag, pivots, work, allowzeropivot, &zeropivotdetected));
6897b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
6904c16a6a6SHong Zhang 
6919371c9d4SSatish Balay     jmin = bi[k];
6929371c9d4SSatish Balay     jmax = bi[k + 1];
6934c16a6a6SHong Zhang     if (jmin < jmax) {
6944c16a6a6SHong Zhang       for (j = jmin; j < jmax; j++) {
6954c16a6a6SHong Zhang         vj       = bj[j]; /* block col. index of U */
6964c16a6a6SHong Zhang         u        = ba + j * bs2;
6974c16a6a6SHong Zhang         rtmp_ptr = rtmp + vj * bs2;
6984c16a6a6SHong Zhang         for (k1 = 0; k1 < bs2; k1++) {
6994c16a6a6SHong Zhang           *u++        = *rtmp_ptr;
7004c16a6a6SHong Zhang           *rtmp_ptr++ = 0.0;
7014c16a6a6SHong Zhang         }
7024c16a6a6SHong Zhang       }
7034c16a6a6SHong Zhang 
7044c16a6a6SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
7054c16a6a6SHong Zhang       il[k] = jmin;
7064c16a6a6SHong Zhang       i     = bj[jmin];
7079371c9d4SSatish Balay       jl[k] = jl[i];
7089371c9d4SSatish Balay       jl[i] = k;
7094c16a6a6SHong Zhang     }
7104c16a6a6SHong Zhang   }
7114c16a6a6SHong Zhang 
7129566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
7139566063dSJacob Faibussowitsch   PetscCall(PetscFree2(il, jl));
7149566063dSJacob Faibussowitsch   PetscCall(PetscFree3(dk, uik, work));
7159566063dSJacob Faibussowitsch   PetscCall(PetscFree(pivots));
7161baa6e33SBarry Smith   if (a->permute) PetscCall(PetscFree(aa));
7174c16a6a6SHong Zhang 
7189566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(perm, &perm_ptr));
71926fbe8dcSKarl Rupp 
7204f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_N_inplace;
7214f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_inplace;
7224f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_inplace;
7234f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_inplace;
724db4efbfdSBarry Smith 
7254c16a6a6SHong Zhang   C->assembled    = PETSC_TRUE;
7264c16a6a6SHong Zhang   C->preallocated = PETSC_TRUE;
72726fbe8dcSKarl Rupp 
7289566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.3333 * bs * bs2 * b->mbs)); /* from inverting diagonal blocks */
7293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7304c16a6a6SHong Zhang }
731d4edadd4SHong Zhang 
732d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat C, Mat A, const MatFactorInfo *info)
733d71ae5a4SJacob Faibussowitsch {
734671cb588SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data, *b = (Mat_SeqSBAIJ *)C->data;
73513f74950SBarry Smith   PetscInt      i, j, mbs = a->mbs, *bi = b->i, *bj = b->j;
73613f74950SBarry Smith   PetscInt     *ai, *aj, k, k1, jmin, jmax, *jl, *il, vj, nexti, ili;
7373bc0b13bSBarry Smith   PetscInt      bs = A->rmap->bs, bs2 = a->bs2;
738671cb588SHong Zhang   MatScalar    *ba = b->a, *aa, *ap, *dk, *uik;
739671cb588SHong Zhang   MatScalar    *u, *diag, *rtmp, *rtmp_ptr;
74028de702eSHong Zhang   MatScalar    *work;
74113f74950SBarry Smith   PetscInt     *pivots;
7425f8bbccaSHong Zhang   PetscBool     allowzeropivot, zeropivotdetected;
743671cb588SHong Zhang 
744671cb588SHong Zhang   PetscFunctionBegin;
7459566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(bs2 * mbs, &rtmp));
7469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(mbs, &il, mbs, &jl));
7476df5ee2eSHong Zhang   il[0] = 0;
7486df5ee2eSHong Zhang   for (i = 0; i < mbs; i++) jl[i] = mbs;
7496df5ee2eSHong Zhang 
7509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(bs2, &dk, bs2, &uik, bs, &work));
7519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs, &pivots));
7525f8bbccaSHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
753671cb588SHong Zhang 
7549371c9d4SSatish Balay   ai = a->i;
7559371c9d4SSatish Balay   aj = a->j;
7569371c9d4SSatish Balay   aa = a->a;
757671cb588SHong Zhang 
758671cb588SHong Zhang   /* for each row k */
759671cb588SHong Zhang   for (k = 0; k < mbs; k++) {
760671cb588SHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
7619371c9d4SSatish Balay     jmin = ai[k];
7629371c9d4SSatish Balay     jmax = ai[k + 1];
763671cb588SHong Zhang     ap   = aa + jmin * bs2;
764671cb588SHong Zhang     for (j = jmin; j < jmax; j++) {
765671cb588SHong Zhang       vj       = aj[j]; /* block col. index */
766671cb588SHong Zhang       rtmp_ptr = rtmp + vj * bs2;
767671cb588SHong Zhang       for (i = 0; i < bs2; i++) *rtmp_ptr++ = *ap++;
768671cb588SHong Zhang     }
769671cb588SHong Zhang 
770671cb588SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
7719566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dk, rtmp + k * bs2, bs2));
772671cb588SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
773671cb588SHong Zhang 
774057f5ba7SHong Zhang     while (i < k) {
775671cb588SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
776671cb588SHong Zhang 
777671cb588SHong Zhang       /* compute multiplier */
778671cb588SHong Zhang       ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
779671cb588SHong Zhang 
780671cb588SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
781671cb588SHong Zhang       diag = ba + i * bs2;
782671cb588SHong Zhang       u    = ba + ili * bs2;
7839566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(uik, bs2));
78496b95a6bSBarry Smith       PetscKernel_A_gets_A_minus_B_times_C(bs, uik, diag, u);
785671cb588SHong Zhang 
786671cb588SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
78796b95a6bSBarry Smith       PetscKernel_A_gets_A_plus_Btranspose_times_C(bs, dk, uik, u);
7889566063dSJacob Faibussowitsch       PetscCall(PetscLogFlops(2.0 * bs * bs2));
789671cb588SHong Zhang 
790671cb588SHong Zhang       /* update -U(i,k) */
7919566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(ba + ili * bs2, uik, bs2));
792671cb588SHong Zhang 
793671cb588SHong Zhang       /* add multiple of row i to k-th row ... */
7949371c9d4SSatish Balay       jmin = ili + 1;
7959371c9d4SSatish Balay       jmax = bi[i + 1];
796671cb588SHong Zhang       if (jmin < jmax) {
797671cb588SHong Zhang         for (j = jmin; j < jmax; j++) {
798671cb588SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
799671cb588SHong Zhang           rtmp_ptr = rtmp + bj[j] * bs2;
800671cb588SHong Zhang           u        = ba + j * bs2;
80196b95a6bSBarry Smith           PetscKernel_A_gets_A_plus_Btranspose_times_C(bs, rtmp_ptr, uik, u);
802671cb588SHong Zhang         }
8039566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(2.0 * bs * bs2 * (jmax - jmin)));
804671cb588SHong Zhang 
805671cb588SHong Zhang         /* ... add i to row list for next nonzero entry */
806671cb588SHong Zhang         il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
807671cb588SHong Zhang         j     = bj[jmin];
8089371c9d4SSatish Balay         jl[i] = jl[j];
8099371c9d4SSatish Balay         jl[j] = i; /* update jl */
810671cb588SHong Zhang       }
811671cb588SHong Zhang       i = nexti;
812671cb588SHong Zhang     }
813671cb588SHong Zhang 
814671cb588SHong Zhang     /* save nonzero entries in k-th row of U ... */
815671cb588SHong Zhang 
816671cb588SHong Zhang     /* invert diagonal block */
817671cb588SHong Zhang     diag = ba + k * bs2;
8189566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(diag, dk, bs2));
8195f8bbccaSHong Zhang 
8209566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A(bs, diag, pivots, work, allowzeropivot, &zeropivotdetected));
8217b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
822671cb588SHong Zhang 
8239371c9d4SSatish Balay     jmin = bi[k];
8249371c9d4SSatish Balay     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];
8399371c9d4SSatish Balay       jl[k] = jl[i];
8409371c9d4SSatish Balay       jl[i] = k;
841671cb588SHong Zhang     }
842671cb588SHong Zhang   }
843671cb588SHong Zhang 
8449566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
8459566063dSJacob Faibussowitsch   PetscCall(PetscFree2(il, jl));
8469566063dSJacob Faibussowitsch   PetscCall(PetscFree3(dk, uik, work));
8479566063dSJacob Faibussowitsch   PetscCall(PetscFree(pivots));
848671cb588SHong Zhang 
8494f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8504f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8514f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8524f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
853671cb588SHong Zhang   C->assembled           = PETSC_TRUE;
854671cb588SHong Zhang   C->preallocated        = PETSC_TRUE;
85526fbe8dcSKarl Rupp 
8569566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.3333 * bs * bs2 * b->mbs)); /* from inverting diagonal blocks */
8573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
858671cb588SHong Zhang }
859671cb588SHong Zhang 
86049b5e25fSSatish Balay /*
861fcf159c0SHong Zhang     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
862cc0c071aSHong Zhang     Version for blocks 2 by 2.
86349b5e25fSSatish Balay */
864d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat C, Mat A, const MatFactorInfo *info)
865d71ae5a4SJacob Faibussowitsch {
86691602c66SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ *)A->data, *b = (Mat_SeqSBAIJ *)C->data;
867cc0c071aSHong Zhang   IS              perm = b->row;
8685d0c19d7SBarry Smith   const PetscInt *ai, *aj, *perm_ptr;
8695d0c19d7SBarry Smith   PetscInt        i, j, mbs = a->mbs, *bi = b->i, *bj = b->j;
8705d0c19d7SBarry Smith   PetscInt       *a2anew, k, k1, jmin, jmax, *jl, *il, vj, nexti, ili;
871d8c74875SBarry Smith   MatScalar      *ba = b->a, *aa, *ap;
872d8c74875SBarry Smith   MatScalar      *u, *diag, *rtmp, *rtmp_ptr, dk[4], uik[4];
87314d2772eSHong Zhang   PetscReal       shift = info->shiftamount;
874a455e926SHong Zhang   PetscBool       allowzeropivot, zeropivotdetected;
87549b5e25fSSatish Balay 
87649b5e25fSSatish Balay   PetscFunctionBegin;
8770164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
8780164db54SHong Zhang 
87991602c66SHong Zhang   /* initialization */
88091602c66SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
88191602c66SHong Zhang      window U(0:k, k:mbs-1).
88291602c66SHong Zhang      jl:    list of rows to be added to uneliminated rows
88391602c66SHong Zhang             i>= k: jl(i) is the first row to be added to row i
88491602c66SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
88591602c66SHong Zhang             jl(i) = mbs indicates the end of a list
88691602c66SHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
88791602c66SHong Zhang             row i of U */
8889566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(4 * mbs, &rtmp));
8899566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(mbs, &il, mbs, &jl));
8906df5ee2eSHong Zhang   il[0] = 0;
8916df5ee2eSHong Zhang   for (i = 0; i < mbs; i++) jl[i] = mbs;
8926df5ee2eSHong Zhang 
8939566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(perm, &perm_ptr));
894cc0c071aSHong Zhang 
895cc0c071aSHong Zhang   /* check permutation */
896cc0c071aSHong Zhang   if (!a->permute) {
8979371c9d4SSatish Balay     ai = a->i;
8989371c9d4SSatish Balay     aj = a->j;
8999371c9d4SSatish Balay     aa = a->a;
900cc0c071aSHong Zhang   } else {
9019371c9d4SSatish Balay     ai = a->inew;
9029371c9d4SSatish Balay     aj = a->jnew;
9039566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(4 * ai[mbs], &aa));
9049566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(aa, a->a, 4 * ai[mbs]));
9059566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(ai[mbs], &a2anew));
9069566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(a2anew, a->a2anew, ai[mbs]));
907cc0c071aSHong Zhang 
908cc0c071aSHong Zhang     for (i = 0; i < mbs; i++) {
9099371c9d4SSatish Balay       jmin = ai[i];
9109371c9d4SSatish Balay       jmax = ai[i + 1];
911cc0c071aSHong Zhang       for (j = jmin; j < jmax; j++) {
912cc0c071aSHong Zhang         while (a2anew[j] != j) {
9139371c9d4SSatish Balay           k         = a2anew[j];
9149371c9d4SSatish Balay           a2anew[j] = a2anew[k];
9159371c9d4SSatish Balay           a2anew[k] = k;
916cc0c071aSHong Zhang           for (k1 = 0; k1 < 4; k1++) {
917cc0c071aSHong Zhang             dk[k1]         = aa[k * 4 + k1];
918cc0c071aSHong Zhang             aa[k * 4 + k1] = aa[j * 4 + k1];
919cc0c071aSHong Zhang             aa[j * 4 + k1] = dk[k1];
920cc0c071aSHong Zhang           }
921cc0c071aSHong Zhang         }
9225e116b59SBarry Smith         /* transform column-oriented blocks that lie in the lower triangle to row-oriented blocks */
923cc0c071aSHong Zhang         if (i > aj[j]) {
924cc0c071aSHong Zhang           ap    = aa + j * 4; /* ptr to the beginning of the block */
925cc0c071aSHong Zhang           dk[1] = ap[1];      /* swap ap[1] and ap[2] */
926cc0c071aSHong Zhang           ap[1] = ap[2];
927cc0c071aSHong Zhang           ap[2] = dk[1];
928cc0c071aSHong Zhang         }
929cc0c071aSHong Zhang       }
930cc0c071aSHong Zhang     }
9319566063dSJacob Faibussowitsch     PetscCall(PetscFree(a2anew));
932cc0c071aSHong Zhang   }
9333845f261SHong Zhang 
93491602c66SHong Zhang   /* for each row k */
93591602c66SHong Zhang   for (k = 0; k < mbs; k++) {
93691602c66SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
9379371c9d4SSatish Balay     jmin = ai[perm_ptr[k]];
9389371c9d4SSatish Balay     jmax = ai[perm_ptr[k] + 1];
939cc0c071aSHong Zhang     ap   = aa + jmin * 4;
94091602c66SHong Zhang     for (j = jmin; j < jmax; j++) {
941cc0c071aSHong Zhang       vj       = perm_ptr[aj[j]]; /* block col. index */
942cc0c071aSHong Zhang       rtmp_ptr = rtmp + vj * 4;
943cc0c071aSHong Zhang       for (i = 0; i < 4; i++) *rtmp_ptr++ = *ap++;
94491602c66SHong Zhang     }
94591602c66SHong Zhang 
94691602c66SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
9479566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dk, rtmp + k * 4, 4));
94891602c66SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
94991602c66SHong Zhang 
950057f5ba7SHong Zhang     while (i < k) {
95191602c66SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
95291602c66SHong Zhang 
9533845f261SHong Zhang       /* compute multiplier */
95491602c66SHong Zhang       ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
9553845f261SHong Zhang 
9563845f261SHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
957cc0c071aSHong Zhang       diag   = ba + i * 4;
958cc0c071aSHong Zhang       u      = ba + ili * 4;
959cc0c071aSHong Zhang       uik[0] = -(diag[0] * u[0] + diag[2] * u[1]);
960cc0c071aSHong Zhang       uik[1] = -(diag[1] * u[0] + diag[3] * u[1]);
961cc0c071aSHong Zhang       uik[2] = -(diag[0] * u[2] + diag[2] * u[3]);
962cc0c071aSHong Zhang       uik[3] = -(diag[1] * u[2] + diag[3] * u[3]);
9633845f261SHong Zhang 
9643845f261SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
965cc0c071aSHong Zhang       dk[0] += uik[0] * u[0] + uik[1] * u[1];
966cc0c071aSHong Zhang       dk[1] += uik[2] * u[0] + uik[3] * u[1];
967cc0c071aSHong Zhang       dk[2] += uik[0] * u[2] + uik[1] * u[3];
968cc0c071aSHong Zhang       dk[3] += uik[2] * u[2] + uik[3] * u[3];
9693845f261SHong Zhang 
9709566063dSJacob Faibussowitsch       PetscCall(PetscLogFlops(16.0 * 2.0));
971187a9f4bSHong Zhang 
9723845f261SHong Zhang       /* update -U(i,k): ba[ili] = uik */
9739566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(ba + ili * 4, uik, 4));
97491602c66SHong Zhang 
97591602c66SHong Zhang       /* add multiple of row i to k-th row ... */
9769371c9d4SSatish Balay       jmin = ili + 1;
9779371c9d4SSatish Balay       jmax = bi[i + 1];
97891602c66SHong Zhang       if (jmin < jmax) {
9793845f261SHong Zhang         for (j = jmin; j < jmax; j++) {
9803845f261SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
981cc0c071aSHong Zhang           rtmp_ptr = rtmp + bj[j] * 4;
982cc0c071aSHong Zhang           u        = ba + j * 4;
983cc0c071aSHong Zhang           rtmp_ptr[0] += uik[0] * u[0] + uik[1] * u[1];
984cc0c071aSHong Zhang           rtmp_ptr[1] += uik[2] * u[0] + uik[3] * u[1];
985cc0c071aSHong Zhang           rtmp_ptr[2] += uik[0] * u[2] + uik[1] * u[3];
986cc0c071aSHong Zhang           rtmp_ptr[3] += uik[2] * u[2] + uik[3] * u[3];
9873845f261SHong Zhang         }
9889566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(16.0 * (jmax - jmin)));
9893845f261SHong Zhang 
99091602c66SHong Zhang         /* ... add i to row list for next nonzero entry */
99191602c66SHong Zhang         il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
99291602c66SHong Zhang         j     = bj[jmin];
9939371c9d4SSatish Balay         jl[i] = jl[j];
9949371c9d4SSatish Balay         jl[j] = i; /* update jl */
99591602c66SHong Zhang       }
996a1723e09SHong Zhang       i = nexti;
99791602c66SHong Zhang     }
998cc0c071aSHong Zhang 
99991602c66SHong Zhang     /* save nonzero entries in k-th row of U ... */
10003845f261SHong Zhang 
1001cc0c071aSHong Zhang     /* invert diagonal block */
1002cc0c071aSHong Zhang     diag = ba + k * 4;
10039566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(diag, dk, 4));
10049566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_2(diag, shift, allowzeropivot, &zeropivotdetected));
10057b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
10063845f261SHong Zhang 
10079371c9d4SSatish Balay     jmin = bi[k];
10089371c9d4SSatish Balay     jmax = bi[k + 1];
100991602c66SHong Zhang     if (jmin < jmax) {
101091602c66SHong Zhang       for (j = jmin; j < jmax; j++) {
1011cc0c071aSHong Zhang         vj       = bj[j]; /* block col. index of U */
1012cc0c071aSHong Zhang         u        = ba + j * 4;
1013cc0c071aSHong Zhang         rtmp_ptr = rtmp + vj * 4;
1014cc0c071aSHong Zhang         for (k1 = 0; k1 < 4; k1++) {
1015cc0c071aSHong Zhang           *u++        = *rtmp_ptr;
1016cc0c071aSHong Zhang           *rtmp_ptr++ = 0.0;
1017cc0c071aSHong Zhang         }
1018cc0c071aSHong Zhang       }
10193845f261SHong Zhang 
102091602c66SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
102191602c66SHong Zhang       il[k] = jmin;
102291602c66SHong Zhang       i     = bj[jmin];
10239371c9d4SSatish Balay       jl[k] = jl[i];
10249371c9d4SSatish Balay       jl[i] = k;
102591602c66SHong Zhang     }
102691602c66SHong Zhang   }
10273845f261SHong Zhang 
10289566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
10299566063dSJacob Faibussowitsch   PetscCall(PetscFree2(il, jl));
10301baa6e33SBarry Smith   if (a->permute) PetscCall(PetscFree(aa));
10319566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(perm, &perm_ptr));
103226fbe8dcSKarl Rupp 
10334f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_2_inplace;
10344f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_inplace;
103549b5e25fSSatish Balay   C->assembled           = PETSC_TRUE;
10365c0bcdfcSHong Zhang   C->preallocated        = PETSC_TRUE;
103726fbe8dcSKarl Rupp 
10389566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.3333 * 8 * b->mbs)); /* from inverting diagonal blocks */
10393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
104049b5e25fSSatish Balay }
104191602c66SHong Zhang 
104249b5e25fSSatish Balay /*
104349b5e25fSSatish Balay       Version for when blocks are 2 by 2 Using natural ordering
104449b5e25fSSatish Balay */
1045d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat C, Mat A, const MatFactorInfo *info)
1046d71ae5a4SJacob Faibussowitsch {
1047ab27746eSHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data, *b = (Mat_SeqSBAIJ *)C->data;
104813f74950SBarry Smith   PetscInt      i, j, mbs = a->mbs, *bi = b->i, *bj = b->j;
104913f74950SBarry Smith   PetscInt     *ai, *aj, k, k1, jmin, jmax, *jl, *il, vj, nexti, ili;
1050d8c74875SBarry Smith   MatScalar    *ba = b->a, *aa, *ap, dk[8], uik[8];
1051ab27746eSHong Zhang   MatScalar    *u, *diag, *rtmp, *rtmp_ptr;
105214d2772eSHong Zhang   PetscReal     shift = info->shiftamount;
1053a455e926SHong Zhang   PetscBool     allowzeropivot, zeropivotdetected;
105449b5e25fSSatish Balay 
105549b5e25fSSatish Balay   PetscFunctionBegin;
10560164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
10570164db54SHong Zhang 
1058ab27746eSHong Zhang   /* initialization */
1059ab27746eSHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
1060ab27746eSHong Zhang      window U(0:k, k:mbs-1).
1061ab27746eSHong Zhang      jl:    list of rows to be added to uneliminated rows
1062ab27746eSHong Zhang             i>= k: jl(i) is the first row to be added to row i
1063ab27746eSHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
1064ab27746eSHong Zhang             jl(i) = mbs indicates the end of a list
1065ab27746eSHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
1066ab27746eSHong Zhang             row i of U */
10679566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(4 * mbs, &rtmp));
10689566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(mbs, &il, mbs, &jl));
10696df5ee2eSHong Zhang   il[0] = 0;
10706df5ee2eSHong Zhang   for (i = 0; i < mbs; i++) jl[i] = mbs;
10716df5ee2eSHong Zhang 
10729371c9d4SSatish Balay   ai = a->i;
10739371c9d4SSatish Balay   aj = a->j;
10749371c9d4SSatish Balay   aa = a->a;
1075ab27746eSHong Zhang 
1076ab27746eSHong Zhang   /* for each row k */
1077ab27746eSHong Zhang   for (k = 0; k < mbs; k++) {
1078ab27746eSHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
10799371c9d4SSatish Balay     jmin = ai[k];
10809371c9d4SSatish Balay     jmax = ai[k + 1];
1081ab27746eSHong Zhang     ap   = aa + jmin * 4;
1082ab27746eSHong Zhang     for (j = jmin; j < jmax; j++) {
1083ab27746eSHong Zhang       vj       = aj[j]; /* block col. index */
1084ab27746eSHong Zhang       rtmp_ptr = rtmp + vj * 4;
1085ab27746eSHong Zhang       for (i = 0; i < 4; i++) *rtmp_ptr++ = *ap++;
108649b5e25fSSatish Balay     }
1087ab27746eSHong Zhang 
1088ab27746eSHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
10899566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dk, rtmp + k * 4, 4));
1090ab27746eSHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
1091ab27746eSHong Zhang 
1092057f5ba7SHong Zhang     while (i < k) {
1093ab27746eSHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
1094ab27746eSHong Zhang 
1095ab27746eSHong Zhang       /* compute multiplier */
1096ab27746eSHong Zhang       ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1097ab27746eSHong Zhang 
1098ab27746eSHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
1099ab27746eSHong Zhang       diag   = ba + i * 4;
1100ab27746eSHong Zhang       u      = ba + ili * 4;
1101ab27746eSHong Zhang       uik[0] = -(diag[0] * u[0] + diag[2] * u[1]);
1102ab27746eSHong Zhang       uik[1] = -(diag[1] * u[0] + diag[3] * u[1]);
1103ab27746eSHong Zhang       uik[2] = -(diag[0] * u[2] + diag[2] * u[3]);
1104ab27746eSHong Zhang       uik[3] = -(diag[1] * u[2] + diag[3] * u[3]);
1105ab27746eSHong Zhang 
1106ab27746eSHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
1107ab27746eSHong Zhang       dk[0] += uik[0] * u[0] + uik[1] * u[1];
1108ab27746eSHong Zhang       dk[1] += uik[2] * u[0] + uik[3] * u[1];
1109ab27746eSHong Zhang       dk[2] += uik[0] * u[2] + uik[1] * u[3];
1110ab27746eSHong Zhang       dk[3] += uik[2] * u[2] + uik[3] * u[3];
1111ab27746eSHong Zhang 
11129566063dSJacob Faibussowitsch       PetscCall(PetscLogFlops(16.0 * 2.0));
1113187a9f4bSHong Zhang 
1114ab27746eSHong Zhang       /* update -U(i,k): ba[ili] = uik */
11159566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(ba + ili * 4, uik, 4));
1116ab27746eSHong Zhang 
1117ab27746eSHong Zhang       /* add multiple of row i to k-th row ... */
11189371c9d4SSatish Balay       jmin = ili + 1;
11199371c9d4SSatish Balay       jmax = bi[i + 1];
1120ab27746eSHong Zhang       if (jmin < jmax) {
1121ab27746eSHong Zhang         for (j = jmin; j < jmax; j++) {
1122ab27746eSHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
1123ab27746eSHong Zhang           rtmp_ptr = rtmp + bj[j] * 4;
1124ab27746eSHong Zhang           u        = ba + j * 4;
1125ab27746eSHong Zhang           rtmp_ptr[0] += uik[0] * u[0] + uik[1] * u[1];
1126ab27746eSHong Zhang           rtmp_ptr[1] += uik[2] * u[0] + uik[3] * u[1];
1127ab27746eSHong Zhang           rtmp_ptr[2] += uik[0] * u[2] + uik[1] * u[3];
1128ab27746eSHong Zhang           rtmp_ptr[3] += uik[2] * u[2] + uik[3] * u[3];
112949b5e25fSSatish Balay         }
11309566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(16.0 * (jmax - jmin)));
1131ab27746eSHong Zhang 
1132ab27746eSHong Zhang         /* ... add i to row list for next nonzero entry */
1133ab27746eSHong Zhang         il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
1134ab27746eSHong Zhang         j     = bj[jmin];
11359371c9d4SSatish Balay         jl[i] = jl[j];
11369371c9d4SSatish Balay         jl[j] = i; /* update jl */
113749b5e25fSSatish Balay       }
1138ab27746eSHong Zhang       i = nexti;
113949b5e25fSSatish Balay     }
1140ab27746eSHong Zhang 
1141ab27746eSHong Zhang     /* save nonzero entries in k-th row of U ... */
1142ab27746eSHong Zhang 
114349b5e25fSSatish Balay     /* invert diagonal block */
1144ab27746eSHong Zhang     diag = ba + k * 4;
11459566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(diag, dk, 4));
11469566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_2(diag, shift, allowzeropivot, &zeropivotdetected));
11477b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1148ab27746eSHong Zhang 
11499371c9d4SSatish Balay     jmin = bi[k];
11509371c9d4SSatish Balay     jmax = bi[k + 1];
1151ab27746eSHong Zhang     if (jmin < jmax) {
1152ab27746eSHong Zhang       for (j = jmin; j < jmax; j++) {
1153ab27746eSHong Zhang         vj       = bj[j]; /* block col. index of U */
1154ab27746eSHong Zhang         u        = ba + j * 4;
1155ab27746eSHong Zhang         rtmp_ptr = rtmp + vj * 4;
1156ab27746eSHong Zhang         for (k1 = 0; k1 < 4; k1++) {
1157ab27746eSHong Zhang           *u++        = *rtmp_ptr;
1158ab27746eSHong Zhang           *rtmp_ptr++ = 0.0;
1159ab27746eSHong Zhang         }
1160ab27746eSHong Zhang       }
1161ab27746eSHong Zhang 
1162ab27746eSHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
1163ab27746eSHong Zhang       il[k] = jmin;
1164ab27746eSHong Zhang       i     = bj[jmin];
11659371c9d4SSatish Balay       jl[k] = jl[i];
11669371c9d4SSatish Balay       jl[i] = k;
1167ab27746eSHong Zhang     }
116849b5e25fSSatish Balay   }
116949b5e25fSSatish Balay 
11709566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
11719566063dSJacob Faibussowitsch   PetscCall(PetscFree2(il, jl));
1172ab27746eSHong Zhang 
11734f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11744f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11754f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11764f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
117749b5e25fSSatish Balay   C->assembled           = PETSC_TRUE;
11785c0bcdfcSHong Zhang   C->preallocated        = PETSC_TRUE;
117926fbe8dcSKarl Rupp 
11809566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.3333 * 8 * b->mbs)); /* from inverting diagonal blocks */
11813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
118249b5e25fSSatish Balay }
118349b5e25fSSatish Balay 
118449b5e25fSSatish Balay /*
118598a8e78dSHong Zhang     Numeric U^T*D*U factorization for SBAIJ format.
118691602c66SHong Zhang     Version for blocks are 1 by 1.
118749b5e25fSSatish Balay */
1188d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace(Mat C, Mat A, const MatFactorInfo *info)
1189d71ae5a4SJacob Faibussowitsch {
119049b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ *)A->data, *b = (Mat_SeqSBAIJ *)C->data;
119149b5e25fSSatish Balay   IS              ip = b->row;
11925d0c19d7SBarry Smith   const PetscInt *ai, *aj, *rip;
11935d0c19d7SBarry Smith   PetscInt       *a2anew, i, j, mbs = a->mbs, *bi = b->i, *bj = b->j, *bcol;
1194997a0949SHong Zhang   PetscInt        k, jmin, jmax, *jl, *il, col, nexti, ili, nz;
1195997a0949SHong Zhang   MatScalar      *rtmp, *ba = b->a, *bval, *aa, dk, uikdi;
11960e95ead3SHong Zhang   PetscReal       rs;
11970e95ead3SHong Zhang   FactorShiftCtx  sctx;
119849b5e25fSSatish Balay 
119949b5e25fSSatish Balay   PetscFunctionBegin;
12000e95ead3SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
12019566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&sctx, sizeof(FactorShiftCtx)));
12023cea4cbeSHong Zhang 
12039566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(ip, &rip));
1204cb718733SHong Zhang   if (!a->permute) {
12059371c9d4SSatish Balay     ai = a->i;
12069371c9d4SSatish Balay     aj = a->j;
12079371c9d4SSatish Balay     aa = a->a;
12082d007305SHong Zhang   } else {
12099371c9d4SSatish Balay     ai = a->inew;
12109371c9d4SSatish Balay     aj = a->jnew;
1211fff829cfSHong Zhang     nz = ai[mbs];
12129566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &aa));
1213fff829cfSHong Zhang     a2anew = a->a2anew;
1214997a0949SHong Zhang     bval   = a->a;
1215ad540459SPierre Jolivet     for (j = 0; j < nz; j++) aa[a2anew[j]] = *(bval++);
12162d007305SHong Zhang   }
12172d007305SHong Zhang 
121891602c66SHong Zhang   /* initialization */
121949b5e25fSSatish Balay   /* il and jl record the first nonzero element in each row of the accessing
122049b5e25fSSatish Balay      window U(0:k, k:mbs-1).
122149b5e25fSSatish Balay      jl:    list of rows to be added to uneliminated rows
122249b5e25fSSatish Balay             i>= k: jl(i) is the first row to be added to row i
122349b5e25fSSatish Balay             i<  k: jl(i) is the row following row i in some list of rows
122449b5e25fSSatish Balay             jl(i) = mbs indicates the end of a list
122549b5e25fSSatish Balay      il(i): points to the first nonzero element in columns k,...,mbs-1 of
122649b5e25fSSatish Balay             row i of U */
12279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(mbs, &rtmp, mbs, &il, mbs, &jl));
1228997a0949SHong Zhang 
1229997a0949SHong Zhang   do {
123007b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
12316df5ee2eSHong Zhang     il[0]         = 0;
123249b5e25fSSatish Balay     for (i = 0; i < mbs; i++) {
12339371c9d4SSatish Balay       rtmp[i] = 0.0;
12349371c9d4SSatish Balay       jl[i]   = mbs;
123549b5e25fSSatish Balay     }
1236997a0949SHong Zhang 
123749b5e25fSSatish Balay     for (k = 0; k < mbs; k++) {
1238997a0949SHong Zhang       /*initialize k-th row by the perm[k]-th row of A */
12399371c9d4SSatish Balay       jmin = ai[rip[k]];
12409371c9d4SSatish Balay       jmax = ai[rip[k] + 1];
12417701de02SHong Zhang       bval = ba + bi[k];
124249b5e25fSSatish Balay       for (j = jmin; j < jmax; j++) {
1243997a0949SHong Zhang         col       = rip[aj[j]];
1244997a0949SHong Zhang         rtmp[col] = aa[j];
12457701de02SHong Zhang         *bval++   = 0.0; /* for in-place factorization */
124649b5e25fSSatish Balay       }
12473cea4cbeSHong Zhang 
12483cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
12493cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
125049b5e25fSSatish Balay 
125191602c66SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
125249b5e25fSSatish Balay       dk = rtmp[k];
125349b5e25fSSatish Balay       i  = jl[k]; /* first row to be added to k_th row  */
125449b5e25fSSatish Balay 
1255057f5ba7SHong Zhang       while (i < k) {
125649b5e25fSSatish Balay         nexti = jl[i]; /* next row to be added to k_th row */
125749b5e25fSSatish Balay 
1258fff829cfSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
125949b5e25fSSatish Balay         ili   = il[i];                /* index of first nonzero element in U(i,k:bms-1) */
1260997a0949SHong Zhang         uikdi = -ba[ili] * ba[bi[i]]; /* diagonal(k) */
126149b5e25fSSatish Balay         dk += uikdi * ba[ili];
1262658e7b3eSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
126349b5e25fSSatish Balay 
1264997a0949SHong Zhang         /* add multiple of row i to k-th row */
12659371c9d4SSatish Balay         jmin = ili + 1;
12669371c9d4SSatish Balay         jmax = bi[i + 1];
126749b5e25fSSatish Balay         if (jmin < jmax) {
126849b5e25fSSatish Balay           for (j = jmin; j < jmax; j++) rtmp[bj[j]] += uikdi * ba[j];
12699566063dSJacob Faibussowitsch           PetscCall(PetscLogFlops(2.0 * (jmax - jmin)));
1270187a9f4bSHong Zhang 
1271fff829cfSHong Zhang           /* update il and jl for row i */
1272fff829cfSHong Zhang           il[i] = jmin;
12739371c9d4SSatish Balay           j     = bj[jmin];
12749371c9d4SSatish Balay           jl[i] = jl[j];
12759371c9d4SSatish Balay           jl[j] = i;
127649b5e25fSSatish Balay         }
1277ab27746eSHong Zhang         i = nexti;
127849b5e25fSSatish Balay       }
127949b5e25fSSatish Balay 
12803cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
12813cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
12823cea4cbeSHong Zhang       rs   = 0.0;
1283997a0949SHong Zhang       jmin = bi[k] + 1;
1284997a0949SHong Zhang       nz   = bi[k + 1] - jmin;
1285997a0949SHong Zhang       if (nz) {
1286997a0949SHong Zhang         bcol = bj + jmin;
1287997a0949SHong Zhang         while (nz--) {
128889f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
128989f845c8SHong Zhang           bcol++;
1290997a0949SHong Zhang         }
1291997a0949SHong Zhang       }
12923cea4cbeSHong Zhang 
12933cea4cbeSHong Zhang       sctx.rs = rs;
12943cea4cbeSHong Zhang       sctx.pv = dk;
12959566063dSJacob Faibussowitsch       PetscCall(MatPivotCheck(C, A, info, &sctx, k));
129607b50cabSHong Zhang       if (sctx.newshift) break; /* sctx.shift_amount is updated */
12970e95ead3SHong Zhang       dk = sctx.pv;
129849b5e25fSSatish Balay 
1299997a0949SHong Zhang       /* copy data into U(k,:) */
130098a8e78dSHong Zhang       ba[bi[k]] = 1.0 / dk; /* U(k,k) */
13019371c9d4SSatish Balay       jmin      = bi[k] + 1;
13029371c9d4SSatish Balay       jmax      = bi[k + 1];
130349b5e25fSSatish Balay       if (jmin < jmax) {
130449b5e25fSSatish Balay         for (j = jmin; j < jmax; j++) {
13059371c9d4SSatish Balay           col       = bj[j];
13069371c9d4SSatish Balay           ba[j]     = rtmp[col];
13079371c9d4SSatish Balay           rtmp[col] = 0.0;
130849b5e25fSSatish Balay         }
130998a8e78dSHong Zhang         /* add the k-th row into il and jl */
131049b5e25fSSatish Balay         il[k] = jmin;
13119371c9d4SSatish Balay         i     = bj[jmin];
13129371c9d4SSatish Balay         jl[k] = jl[i];
13139371c9d4SSatish Balay         jl[i] = k;
131449b5e25fSSatish Balay       }
131549b5e25fSSatish Balay     }
131607b50cabSHong Zhang   } while (sctx.newshift);
13179566063dSJacob Faibussowitsch   PetscCall(PetscFree3(rtmp, il, jl));
13189566063dSJacob Faibussowitsch   if (a->permute) PetscCall(PetscFree(aa));
131949b5e25fSSatish Balay 
13209566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(ip, &rip));
132126fbe8dcSKarl Rupp 
13220a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_inplace;
13234f79d315SHong Zhang   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
13240a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace;
13250a3c351aSHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_inplace;
13260a3c351aSHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_inplace;
132749b5e25fSSatish Balay   C->assembled           = PETSC_TRUE;
13285c0bcdfcSHong Zhang   C->preallocated        = PETSC_TRUE;
132926fbe8dcSKarl Rupp 
13309566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(C->rmap->N));
13313cea4cbeSHong Zhang   if (sctx.nshift) {
1332f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
13339566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shiftnz tries %" PetscInt_FMT ", shift_amount %g\n", sctx.nshift, (double)sctx.shift_amount));
1334f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
13359566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shiftpd tries %" PetscInt_FMT ", shift_amount %g\n", sctx.nshift, (double)sctx.shift_amount));
1336997a0949SHong Zhang     }
1337997a0949SHong Zhang   }
13383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
133949b5e25fSSatish Balay }
134049b5e25fSSatish Balay 
1341671cb588SHong Zhang /*
134280d3d5a7SHong Zhang   Version for when blocks are 1 by 1 Using natural ordering under new datastructure
134380d3d5a7SHong Zhang   Modified from MatCholeskyFactorNumeric_SeqAIJ()
1344671cb588SHong Zhang */
1345d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat B, Mat A, const MatFactorInfo *info)
1346d71ae5a4SJacob Faibussowitsch {
1347d595f711SHong Zhang   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ *)A->data;
13487b056e98SHong Zhang   Mat_SeqSBAIJ  *b = (Mat_SeqSBAIJ *)B->data;
1349d595f711SHong Zhang   PetscInt       i, j, mbs = A->rmap->n, *bi = b->i, *bj = b->j, *bdiag = b->diag, *bjtmp;
1350545dd064SHong Zhang   PetscInt      *ai = a->i, *aj = a->j, *ajtmp;
1351d595f711SHong Zhang   PetscInt       k, jmin, jmax, *c2r, *il, col, nexti, ili, nz;
1352d595f711SHong Zhang   MatScalar     *rtmp, *ba = b->a, *bval, *aa = a->a, dk, uikdi;
1353d90ac83dSHong Zhang   FactorShiftCtx sctx;
1354d595f711SHong Zhang   PetscReal      rs;
1355d595f711SHong Zhang   MatScalar      d, *v;
1356d595f711SHong Zhang 
1357d595f711SHong Zhang   PetscFunctionBegin;
13589566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(mbs, &rtmp, mbs, &il, mbs, &c2r));
1359545dd064SHong Zhang 
1360d595f711SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
13619566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&sctx, sizeof(FactorShiftCtx)));
1362d595f711SHong Zhang 
1363f4db908eSBarry Smith   if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
1364d595f711SHong Zhang     sctx.shift_top = info->zeropivot;
136526fbe8dcSKarl Rupp 
13669566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(rtmp, mbs));
136726fbe8dcSKarl Rupp 
1368d595f711SHong Zhang     for (i = 0; i < mbs; i++) {
1369d595f711SHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
1370d595f711SHong Zhang       d = (aa)[a->diag[i]];
1371545dd064SHong Zhang       rtmp[i] += -PetscRealPart(d); /* diagonal entry */
1372545dd064SHong Zhang       ajtmp = aj + ai[i] + 1;       /* exclude diagonal */
1373545dd064SHong Zhang       v     = aa + ai[i] + 1;
1374545dd064SHong Zhang       nz    = ai[i + 1] - ai[i] - 1;
1375545dd064SHong Zhang       for (j = 0; j < nz; j++) {
1376545dd064SHong Zhang         rtmp[i] += PetscAbsScalar(v[j]);
1377545dd064SHong Zhang         rtmp[ajtmp[j]] += PetscAbsScalar(v[j]);
1378545dd064SHong Zhang       }
13790838c725SBarry Smith       if (PetscRealPart(rtmp[i]) > sctx.shift_top) sctx.shift_top = PetscRealPart(rtmp[i]);
1380d595f711SHong Zhang     }
1381d595f711SHong Zhang     sctx.shift_top *= 1.1;
1382d595f711SHong Zhang     sctx.nshift_max = 5;
1383d595f711SHong Zhang     sctx.shift_lo   = 0.;
1384d595f711SHong Zhang     sctx.shift_hi   = 1.;
1385d595f711SHong Zhang   }
1386d595f711SHong Zhang 
1387d595f711SHong Zhang   /* allocate working arrays
1388d595f711SHong Zhang      c2r: linked list, keep track of pivot rows for a given column. c2r[col]: head of the list for a given col
1389d595f711SHong 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
1390d595f711SHong Zhang   */
1391d595f711SHong Zhang   do {
139207b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
1393d595f711SHong Zhang 
1394d595f711SHong Zhang     for (i = 0; i < mbs; i++) c2r[i] = mbs;
13959e95ef84SSatish Balay     if (mbs) il[0] = 0;
1396d595f711SHong Zhang 
1397d595f711SHong Zhang     for (k = 0; k < mbs; k++) {
1398d595f711SHong Zhang       /* zero rtmp */
1399d595f711SHong Zhang       nz    = bi[k + 1] - bi[k];
1400d595f711SHong Zhang       bjtmp = bj + bi[k];
14017b056e98SHong Zhang       for (j = 0; j < nz; j++) rtmp[bjtmp[j]] = 0.0;
1402d595f711SHong Zhang 
1403d595f711SHong Zhang       /* load in initial unfactored row */
1404d595f711SHong Zhang       bval = ba + bi[k];
14059371c9d4SSatish Balay       jmin = ai[k];
14069371c9d4SSatish Balay       jmax = ai[k + 1];
1407d595f711SHong Zhang       for (j = jmin; j < jmax; j++) {
1408d595f711SHong Zhang         col       = aj[j];
1409d595f711SHong Zhang         rtmp[col] = aa[j];
1410d595f711SHong Zhang         *bval++   = 0.0; /* for in-place factorization */
1411d595f711SHong Zhang       }
1412d595f711SHong Zhang       /* shift the diagonal of the matrix: ZeropivotApply() */
1413d595f711SHong Zhang       rtmp[k] += sctx.shift_amount; /* shift the diagonal of the matrix */
1414d595f711SHong Zhang 
1415d595f711SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1416d595f711SHong Zhang       dk = rtmp[k];
1417d595f711SHong Zhang       i  = c2r[k]; /* first row to be added to k_th row  */
1418d595f711SHong Zhang 
1419d595f711SHong Zhang       while (i < k) {
1420d595f711SHong Zhang         nexti = c2r[i]; /* next row to be added to k_th row */
1421d595f711SHong Zhang 
1422d595f711SHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
1423d595f711SHong Zhang         ili   = il[i];                   /* index of first nonzero element in U(i,k:bms-1) */
1424d595f711SHong Zhang         uikdi = -ba[ili] * ba[bdiag[i]]; /* diagonal(k) */
1425d595f711SHong Zhang         dk += uikdi * ba[ili];           /* update diag[k] */
1426d595f711SHong Zhang         ba[ili] = uikdi;                 /* -U(i,k) */
1427d595f711SHong Zhang 
1428d595f711SHong Zhang         /* add multiple of row i to k-th row */
14299371c9d4SSatish Balay         jmin = ili + 1;
14309371c9d4SSatish Balay         jmax = bi[i + 1];
1431d595f711SHong Zhang         if (jmin < jmax) {
1432d595f711SHong Zhang           for (j = jmin; j < jmax; j++) rtmp[bj[j]] += uikdi * ba[j];
1433d595f711SHong Zhang           /* update il and c2r for row i */
1434d595f711SHong Zhang           il[i]  = jmin;
14359371c9d4SSatish Balay           j      = bj[jmin];
14369371c9d4SSatish Balay           c2r[i] = c2r[j];
14379371c9d4SSatish Balay           c2r[j] = i;
1438d595f711SHong Zhang         }
1439d595f711SHong Zhang         i = nexti;
1440d595f711SHong Zhang       }
1441d595f711SHong Zhang 
1442d595f711SHong Zhang       /* copy data into U(k,:) */
1443d595f711SHong Zhang       rs   = 0.0;
14449371c9d4SSatish Balay       jmin = bi[k];
14459371c9d4SSatish Balay       jmax = bi[k + 1] - 1;
1446d595f711SHong Zhang       if (jmin < jmax) {
1447d595f711SHong Zhang         for (j = jmin; j < jmax; j++) {
14489371c9d4SSatish Balay           col   = bj[j];
14499371c9d4SSatish Balay           ba[j] = rtmp[col];
14509371c9d4SSatish Balay           rs += PetscAbsScalar(ba[j]);
1451d595f711SHong Zhang         }
1452d595f711SHong Zhang         /* add the k-th row into il and c2r */
1453d595f711SHong Zhang         il[k]  = jmin;
14549371c9d4SSatish Balay         i      = bj[jmin];
14559371c9d4SSatish Balay         c2r[k] = c2r[i];
14569371c9d4SSatish Balay         c2r[i] = k;
1457d595f711SHong Zhang       }
1458d595f711SHong Zhang 
1459d595f711SHong Zhang       sctx.rs = rs;
1460d595f711SHong Zhang       sctx.pv = dk;
14619566063dSJacob Faibussowitsch       PetscCall(MatPivotCheck(B, A, info, &sctx, k));
146207b50cabSHong Zhang       if (sctx.newshift) break;
1463d595f711SHong Zhang       dk = sctx.pv;
1464d595f711SHong Zhang 
1465d595f711SHong Zhang       ba[bdiag[k]] = 1.0 / dk; /* U(k,k) */
1466d595f711SHong Zhang     }
146707b50cabSHong Zhang   } while (sctx.newshift);
1468d595f711SHong Zhang 
14699566063dSJacob Faibussowitsch   PetscCall(PetscFree3(rtmp, il, c2r));
1470d595f711SHong Zhang 
1471d595f711SHong Zhang   B->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering;
14729dcdb97aSHong Zhang   B->ops->solves         = MatSolves_SeqSBAIJ_1;
1473d595f711SHong Zhang   B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1474910cf402Sprj-   B->ops->matsolve       = MatMatSolve_SeqSBAIJ_1_NaturalOrdering;
147580d3d5a7SHong Zhang   B->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
147680d3d5a7SHong Zhang   B->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;
1477d595f711SHong Zhang 
14787b056e98SHong Zhang   B->assembled    = PETSC_TRUE;
14797b056e98SHong Zhang   B->preallocated = PETSC_TRUE;
148026fbe8dcSKarl Rupp 
14819566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(B->rmap->n));
1482d595f711SHong Zhang 
1483d595f711SHong Zhang   /* MatPivotView() */
1484d595f711SHong Zhang   if (sctx.nshift) {
1485f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
14869566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shift_pd tries %" PetscInt_FMT ", 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));
1487f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
14889566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shift_nz tries %" PetscInt_FMT ", shift_amount %g\n", sctx.nshift, (double)sctx.shift_amount));
1489f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) {
14909566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shift_inblocks applied %" PetscInt_FMT ", each shift_amount %g\n", sctx.nshift, (double)info->shiftamount));
1491d595f711SHong Zhang     }
1492d595f711SHong Zhang   }
14933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1494d595f711SHong Zhang }
1495d595f711SHong Zhang 
1496d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace(Mat C, Mat A, const MatFactorInfo *info)
1497d71ae5a4SJacob Faibussowitsch {
1498671cb588SHong Zhang   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ *)A->data, *b = (Mat_SeqSBAIJ *)C->data;
149913f74950SBarry Smith   PetscInt       i, j, mbs = a->mbs;
150013f74950SBarry Smith   PetscInt      *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j;
15013cea4cbeSHong Zhang   PetscInt       k, jmin, *jl, *il, nexti, ili, *acol, *bcol, nz;
1502653a6975SHong Zhang   MatScalar     *rtmp, *ba = b->a, *aa = a->a, dk, uikdi, *aval, *bval;
150314d2772eSHong Zhang   PetscReal      rs;
15040e95ead3SHong Zhang   FactorShiftCtx sctx;
1505653a6975SHong Zhang 
1506653a6975SHong Zhang   PetscFunctionBegin;
15070e95ead3SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
15089566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&sctx, sizeof(FactorShiftCtx)));
15090e95ead3SHong Zhang 
1510653a6975SHong Zhang   /* initialization */
1511653a6975SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
1512653a6975SHong Zhang      window U(0:k, k:mbs-1).
1513653a6975SHong Zhang      jl:    list of rows to be added to uneliminated rows
1514653a6975SHong Zhang             i>= k: jl(i) is the first row to be added to row i
1515653a6975SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
1516653a6975SHong Zhang             jl(i) = mbs indicates the end of a list
1517653a6975SHong Zhang      il(i): points to the first nonzero element in U(i,k:mbs-1)
1518653a6975SHong Zhang   */
15199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs, &rtmp));
15209566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(mbs, &il, mbs, &jl));
1521b00f7748SHong Zhang 
1522b00f7748SHong Zhang   do {
152307b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
15246df5ee2eSHong Zhang     il[0]         = 0;
1525653a6975SHong Zhang     for (i = 0; i < mbs; i++) {
15269371c9d4SSatish Balay       rtmp[i] = 0.0;
15279371c9d4SSatish Balay       jl[i]   = mbs;
1528653a6975SHong Zhang     }
1529653a6975SHong Zhang 
1530997a0949SHong Zhang     for (k = 0; k < mbs; k++) {
1531653a6975SHong Zhang       /*initialize k-th row with elements nonzero in row perm(k) of A */
1532653a6975SHong Zhang       nz   = ai[k + 1] - ai[k];
1533653a6975SHong Zhang       acol = aj + ai[k];
1534653a6975SHong Zhang       aval = aa + ai[k];
1535653a6975SHong Zhang       bval = ba + bi[k];
1536653a6975SHong Zhang       while (nz--) {
1537653a6975SHong Zhang         rtmp[*acol++] = *aval++;
1538653a6975SHong Zhang         *bval++       = 0.0; /* for in-place factorization */
1539653a6975SHong Zhang       }
15403cea4cbeSHong Zhang 
15413cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
15423cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1543653a6975SHong Zhang 
1544653a6975SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1545653a6975SHong Zhang       dk = rtmp[k];
1546653a6975SHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
1547653a6975SHong Zhang 
1548653a6975SHong Zhang       while (i < k) {
1549653a6975SHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
1550653a6975SHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
1551653a6975SHong Zhang         ili   = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1552653a6975SHong Zhang         uikdi = -ba[ili] * ba[bi[i]];
1553653a6975SHong Zhang         dk += uikdi * ba[ili];
1554653a6975SHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
1555653a6975SHong Zhang 
1556653a6975SHong Zhang         /* add multiple of row i to k-th row ... */
1557653a6975SHong Zhang         jmin = ili + 1;
1558653a6975SHong Zhang         nz   = bi[i + 1] - jmin;
1559653a6975SHong Zhang         if (nz > 0) {
1560653a6975SHong Zhang           bcol = bj + jmin;
1561653a6975SHong Zhang           bval = ba + jmin;
15629566063dSJacob Faibussowitsch           PetscCall(PetscLogFlops(2.0 * nz));
1563653a6975SHong Zhang           while (nz--) rtmp[*bcol++] += uikdi * (*bval++);
1564187a9f4bSHong Zhang 
1565997a0949SHong Zhang           /* update il and jl for i-th row */
1566997a0949SHong Zhang           il[i] = jmin;
15679371c9d4SSatish Balay           j     = bj[jmin];
15689371c9d4SSatish Balay           jl[i] = jl[j];
15699371c9d4SSatish Balay           jl[j] = i;
1570653a6975SHong Zhang         }
1571653a6975SHong Zhang         i = nexti;
1572653a6975SHong Zhang       }
1573653a6975SHong Zhang 
15743cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
15753cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
15763cea4cbeSHong Zhang       rs   = 0.0;
157721c26570Svictorle       jmin = bi[k] + 1;
157821c26570Svictorle       nz   = bi[k + 1] - jmin;
157921c26570Svictorle       if (nz) {
158021c26570Svictorle         bcol = bj + jmin;
158121c26570Svictorle         while (nz--) {
158289f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
158389f845c8SHong Zhang           bcol++;
158421c26570Svictorle         }
158521c26570Svictorle       }
15863cea4cbeSHong Zhang 
15873cea4cbeSHong Zhang       sctx.rs = rs;
15883cea4cbeSHong Zhang       sctx.pv = dk;
15899566063dSJacob Faibussowitsch       PetscCall(MatPivotCheck(C, A, info, &sctx, k));
159007b50cabSHong Zhang       if (sctx.newshift) break; /* sctx.shift_amount is updated */
15910e95ead3SHong Zhang       dk = sctx.pv;
1592653a6975SHong Zhang 
1593997a0949SHong Zhang       /* copy data into U(k,:) */
1594653a6975SHong Zhang       ba[bi[k]] = 1.0 / dk;
1595653a6975SHong Zhang       jmin      = bi[k] + 1;
1596653a6975SHong Zhang       nz        = bi[k + 1] - jmin;
1597653a6975SHong Zhang       if (nz) {
1598653a6975SHong Zhang         bcol = bj + jmin;
1599653a6975SHong Zhang         bval = ba + jmin;
1600653a6975SHong Zhang         while (nz--) {
1601653a6975SHong Zhang           *bval++       = rtmp[*bcol];
1602653a6975SHong Zhang           rtmp[*bcol++] = 0.0;
1603653a6975SHong Zhang         }
1604997a0949SHong Zhang         /* add k-th row into il and jl */
1605653a6975SHong Zhang         il[k] = jmin;
16069371c9d4SSatish Balay         i     = bj[jmin];
16079371c9d4SSatish Balay         jl[k] = jl[i];
16089371c9d4SSatish Balay         jl[i] = k;
1609653a6975SHong Zhang       }
1610b00f7748SHong Zhang     } /* end of for (k = 0; k<mbs; k++) */
161107b50cabSHong Zhang   } while (sctx.newshift);
16129566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
16139566063dSJacob Faibussowitsch   PetscCall(PetscFree2(il, jl));
1614653a6975SHong Zhang 
16150a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
16164f79d315SHong Zhang   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
16170a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
16180a3c351aSHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
16190a3c351aSHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1620db4efbfdSBarry Smith 
1621653a6975SHong Zhang   C->assembled    = PETSC_TRUE;
1622653a6975SHong Zhang   C->preallocated = PETSC_TRUE;
162326fbe8dcSKarl Rupp 
16249566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(C->rmap->N));
16253cea4cbeSHong Zhang   if (sctx.nshift) {
1626f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
16279566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shiftnz tries %" PetscInt_FMT ", shift_amount %g\n", sctx.nshift, (double)sctx.shift_amount));
1628f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
16299566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shiftpd tries %" PetscInt_FMT ", shift_amount %g\n", sctx.nshift, (double)sctx.shift_amount));
1630b00f7748SHong Zhang     }
1631fb10cecfSBarry Smith   }
16323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1633653a6975SHong Zhang }
1634653a6975SHong Zhang 
1635d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A, IS perm, const MatFactorInfo *info)
1636d71ae5a4SJacob Faibussowitsch {
163749b5e25fSSatish Balay   Mat C;
163849b5e25fSSatish Balay 
163949b5e25fSSatish Balay   PetscFunctionBegin;
16409566063dSJacob Faibussowitsch   PetscCall(MatGetFactor(A, "petsc", MAT_FACTOR_CHOLESKY, &C));
16419566063dSJacob Faibussowitsch   PetscCall(MatCholeskyFactorSymbolic(C, A, perm, info));
16429566063dSJacob Faibussowitsch   PetscCall(MatCholeskyFactorNumeric(C, A, info));
164326fbe8dcSKarl Rupp 
1644db4efbfdSBarry Smith   A->ops->solve          = C->ops->solve;
1645db4efbfdSBarry Smith   A->ops->solvetranspose = C->ops->solvetranspose;
164626fbe8dcSKarl Rupp 
16479566063dSJacob Faibussowitsch   PetscCall(MatHeaderMerge(A, &C));
16483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
164949b5e25fSSatish Balay }
1650