xref: /petsc/src/mat/impls/baij/seq/baijfact.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1be1d678aSKris Buschelman 
25d517e7eSBarry Smith /*
3ec3a8b7bSBarry Smith     Factorization code for BAIJ format.
45d517e7eSBarry Smith */
5c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
6af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h>
7ec3a8b7bSBarry Smith 
89371c9d4SSatish Balay PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat B, Mat A, const MatFactorInfo *info) {
9b588c5a2SHong Zhang   Mat             C = B;
10b588c5a2SHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
11b588c5a2SHong Zhang   IS              isrow = b->row, isicol = b->icol;
125a586d82SBarry Smith   const PetscInt *r, *ic;
13bbd65245SShri Abhyankar   PetscInt        i, j, k, nz, nzL, row, *pj;
14bbd65245SShri Abhyankar   const PetscInt  n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j, bs2 = a->bs2;
15bbd65245SShri Abhyankar   const PetscInt *ajtmp, *bjtmp, *bdiag = b->diag;
16bbd65245SShri Abhyankar   MatScalar      *rtmp, *pc, *mwork, *pv;
17bbd65245SShri Abhyankar   MatScalar      *aa = a->a, *v;
18bbd65245SShri Abhyankar   PetscInt        flg;
19182b8fbaSHong Zhang   PetscReal       shift = info->shiftamount;
20a455e926SHong Zhang   PetscBool       allowzeropivot, zeropivotdetected;
21b588c5a2SHong Zhang 
22b588c5a2SHong Zhang   PetscFunctionBegin;
239566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
249566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
250164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
26b588c5a2SHong Zhang 
274c000e2eSHong Zhang   /* generate work space needed by the factorization */
289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork));
299566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(rtmp, bs2 * n));
30b588c5a2SHong Zhang 
31b588c5a2SHong Zhang   for (i = 0; i < n; i++) {
32b588c5a2SHong Zhang     /* zero rtmp */
33b588c5a2SHong Zhang     /* L part */
34b588c5a2SHong Zhang     nz    = bi[i + 1] - bi[i];
35b588c5a2SHong Zhang     bjtmp = bj + bi[i];
36*48a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
37b588c5a2SHong Zhang 
38b588c5a2SHong Zhang     /* U part */
390c4413a7SShri Abhyankar     nz    = bdiag[i] - bdiag[i + 1];
400c4413a7SShri Abhyankar     bjtmp = bj + bdiag[i + 1] + 1;
41*48a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
420c4413a7SShri Abhyankar 
430c4413a7SShri Abhyankar     /* load in initial (unfactored row) */
440c4413a7SShri Abhyankar     nz    = ai[r[i] + 1] - ai[r[i]];
450c4413a7SShri Abhyankar     ajtmp = aj + ai[r[i]];
460c4413a7SShri Abhyankar     v     = aa + bs2 * ai[r[i]];
47*48a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ic[ajtmp[j]], v + bs2 * j, bs2));
480c4413a7SShri Abhyankar 
490c4413a7SShri Abhyankar     /* elimination */
500c4413a7SShri Abhyankar     bjtmp = bj + bi[i];
510c4413a7SShri Abhyankar     nzL   = bi[i + 1] - bi[i];
520c4413a7SShri Abhyankar     for (k = 0; k < nzL; k++) {
530c4413a7SShri Abhyankar       row = bjtmp[k];
540c4413a7SShri Abhyankar       pc  = rtmp + bs2 * row;
55c35f09e5SBarry Smith       for (flg = 0, j = 0; j < bs2; j++) {
56c35f09e5SBarry Smith         if (pc[j] != (PetscScalar)0.0) {
57c35f09e5SBarry Smith           flg = 1;
58c35f09e5SBarry Smith           break;
59c35f09e5SBarry Smith         }
60c35f09e5SBarry Smith       }
610c4413a7SShri Abhyankar       if (flg) {
620c4413a7SShri Abhyankar         pv = b->a + bs2 * bdiag[row];
6396b95a6bSBarry Smith         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
649566063dSJacob Faibussowitsch         PetscCall(PetscKernel_A_gets_A_times_B_2(pc, pv, mwork));
650c4413a7SShri Abhyankar 
66a5b23f4aSJose E. Roman         pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */
670c4413a7SShri Abhyankar         pv = b->a + bs2 * (bdiag[row + 1] + 1);
680c4413a7SShri Abhyankar         nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */
690c4413a7SShri Abhyankar         for (j = 0; j < nz; j++) {
7096b95a6bSBarry Smith           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
710c4413a7SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
720c4413a7SShri Abhyankar           v = rtmp + 4 * pj[j];
739566063dSJacob Faibussowitsch           PetscCall(PetscKernel_A_gets_A_minus_B_times_C_2(v, pc, pv));
740c4413a7SShri Abhyankar           pv += 4;
750c4413a7SShri Abhyankar         }
769566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(16.0 * nz + 12)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
770c4413a7SShri Abhyankar       }
780c4413a7SShri Abhyankar     }
790c4413a7SShri Abhyankar 
800c4413a7SShri Abhyankar     /* finished row so stick it into b->a */
810c4413a7SShri Abhyankar     /* L part */
820c4413a7SShri Abhyankar     pv = b->a + bs2 * bi[i];
830c4413a7SShri Abhyankar     pj = b->j + bi[i];
840c4413a7SShri Abhyankar     nz = bi[i + 1] - bi[i];
85*48a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
860c4413a7SShri Abhyankar 
87a5b23f4aSJose E. Roman     /* Mark diagonal and invert diagonal for simpler triangular solves */
880c4413a7SShri Abhyankar     pv = b->a + bs2 * bdiag[i];
890c4413a7SShri Abhyankar     pj = b->j + bdiag[i];
909566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2));
919566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_2(pv, shift, allowzeropivot, &zeropivotdetected));
927b6c816cSBarry Smith     if (zeropivotdetected) B->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
930c4413a7SShri Abhyankar 
940c4413a7SShri Abhyankar     /* U part */
950c4413a7SShri Abhyankar     pv = b->a + bs2 * (bdiag[i + 1] + 1);
960c4413a7SShri Abhyankar     pj = b->j + bdiag[i + 1] + 1;
970c4413a7SShri Abhyankar     nz = bdiag[i] - bdiag[i + 1] - 1;
98*48a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
990c4413a7SShri Abhyankar   }
1000c4413a7SShri Abhyankar 
1019566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rtmp, mwork));
1029566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
1039566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
10426fbe8dcSKarl Rupp 
1054dd39f65SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_2;
1064dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2;
1070c4413a7SShri Abhyankar   C->assembled           = PETSC_TRUE;
10826fbe8dcSKarl Rupp 
1099566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 2 * 2 * 2 * n)); /* from inverting diagonal blocks */
1100c4413a7SShri Abhyankar   PetscFunctionReturn(0);
1110c4413a7SShri Abhyankar }
1120c4413a7SShri Abhyankar 
1139371c9d4SSatish Balay PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat B, Mat A, const MatFactorInfo *info) {
1144c000e2eSHong Zhang   Mat             C = B;
1154c000e2eSHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
116bbd65245SShri Abhyankar   PetscInt        i, j, k, nz, nzL, row, *pj;
117bbd65245SShri Abhyankar   const PetscInt  n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j, bs2 = a->bs2;
118bbd65245SShri Abhyankar   const PetscInt *ajtmp, *bjtmp, *bdiag = b->diag;
119bbd65245SShri Abhyankar   MatScalar      *rtmp, *pc, *mwork, *pv;
120bbd65245SShri Abhyankar   MatScalar      *aa = a->a, *v;
121bbd65245SShri Abhyankar   PetscInt        flg;
122182b8fbaSHong Zhang   PetscReal       shift = info->shiftamount;
123a455e926SHong Zhang   PetscBool       allowzeropivot, zeropivotdetected;
1244c000e2eSHong Zhang 
1254c000e2eSHong Zhang   PetscFunctionBegin;
1260164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
1270164db54SHong Zhang 
1284c000e2eSHong Zhang   /* generate work space needed by the factorization */
1299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork));
1309566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(rtmp, bs2 * n));
1314c000e2eSHong Zhang 
1324c000e2eSHong Zhang   for (i = 0; i < n; i++) {
1334c000e2eSHong Zhang     /* zero rtmp */
1344c000e2eSHong Zhang     /* L part */
1354c000e2eSHong Zhang     nz    = bi[i + 1] - bi[i];
1364c000e2eSHong Zhang     bjtmp = bj + bi[i];
137*48a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
1384c000e2eSHong Zhang 
1394c000e2eSHong Zhang     /* U part */
140b2b2dd24SShri Abhyankar     nz    = bdiag[i] - bdiag[i + 1];
141b2b2dd24SShri Abhyankar     bjtmp = bj + bdiag[i + 1] + 1;
142*48a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
143b2b2dd24SShri Abhyankar 
144b2b2dd24SShri Abhyankar     /* load in initial (unfactored row) */
145b2b2dd24SShri Abhyankar     nz    = ai[i + 1] - ai[i];
146b2b2dd24SShri Abhyankar     ajtmp = aj + ai[i];
147b2b2dd24SShri Abhyankar     v     = aa + bs2 * ai[i];
148*48a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ajtmp[j], v + bs2 * j, bs2));
149b2b2dd24SShri Abhyankar 
150b2b2dd24SShri Abhyankar     /* elimination */
151b2b2dd24SShri Abhyankar     bjtmp = bj + bi[i];
152b2b2dd24SShri Abhyankar     nzL   = bi[i + 1] - bi[i];
153b2b2dd24SShri Abhyankar     for (k = 0; k < nzL; k++) {
154b2b2dd24SShri Abhyankar       row = bjtmp[k];
155b2b2dd24SShri Abhyankar       pc  = rtmp + bs2 * row;
156c35f09e5SBarry Smith       for (flg = 0, j = 0; j < bs2; j++) {
157c35f09e5SBarry Smith         if (pc[j] != (PetscScalar)0.0) {
158c35f09e5SBarry Smith           flg = 1;
159c35f09e5SBarry Smith           break;
160c35f09e5SBarry Smith         }
161c35f09e5SBarry Smith       }
162b2b2dd24SShri Abhyankar       if (flg) {
163b2b2dd24SShri Abhyankar         pv = b->a + bs2 * bdiag[row];
16496b95a6bSBarry Smith         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
1659566063dSJacob Faibussowitsch         PetscCall(PetscKernel_A_gets_A_times_B_2(pc, pv, mwork));
166b2b2dd24SShri Abhyankar 
167b2b2dd24SShri Abhyankar         pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */
168b2b2dd24SShri Abhyankar         pv = b->a + bs2 * (bdiag[row + 1] + 1);
169b2b2dd24SShri Abhyankar         nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries in U(row,:) excluding diag */
170b2b2dd24SShri Abhyankar         for (j = 0; j < nz; j++) {
17196b95a6bSBarry Smith           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
172b2b2dd24SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
173b2b2dd24SShri Abhyankar           v = rtmp + 4 * pj[j];
1749566063dSJacob Faibussowitsch           PetscCall(PetscKernel_A_gets_A_minus_B_times_C_2(v, pc, pv));
175b2b2dd24SShri Abhyankar           pv += 4;
176b2b2dd24SShri Abhyankar         }
1779566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(16.0 * nz + 12)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
178b2b2dd24SShri Abhyankar       }
179b2b2dd24SShri Abhyankar     }
180b2b2dd24SShri Abhyankar 
181b2b2dd24SShri Abhyankar     /* finished row so stick it into b->a */
182b2b2dd24SShri Abhyankar     /* L part */
183b2b2dd24SShri Abhyankar     pv = b->a + bs2 * bi[i];
184b2b2dd24SShri Abhyankar     pj = b->j + bi[i];
185b2b2dd24SShri Abhyankar     nz = bi[i + 1] - bi[i];
186*48a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
187b2b2dd24SShri Abhyankar 
188a5b23f4aSJose E. Roman     /* Mark diagonal and invert diagonal for simpler triangular solves */
189b2b2dd24SShri Abhyankar     pv = b->a + bs2 * bdiag[i];
190b2b2dd24SShri Abhyankar     pj = b->j + bdiag[i];
1919566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2));
1929566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_2(pv, shift, allowzeropivot, &zeropivotdetected));
1937b6c816cSBarry Smith     if (zeropivotdetected) B->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
194b2b2dd24SShri Abhyankar 
195b2b2dd24SShri Abhyankar     /* U part */
196b2b2dd24SShri Abhyankar     /*
197b2b2dd24SShri Abhyankar     pv = b->a + bs2*bi[2*n-i];
198b2b2dd24SShri Abhyankar     pj = b->j + bi[2*n-i];
199b2b2dd24SShri Abhyankar     nz = bi[2*n-i+1] - bi[2*n-i] - 1;
200b2b2dd24SShri Abhyankar     */
201b2b2dd24SShri Abhyankar     pv = b->a + bs2 * (bdiag[i + 1] + 1);
202b2b2dd24SShri Abhyankar     pj = b->j + bdiag[i + 1] + 1;
203b2b2dd24SShri Abhyankar     nz = bdiag[i] - bdiag[i + 1] - 1;
204*48a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
205b2b2dd24SShri Abhyankar   }
2069566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rtmp, mwork));
207ae3d28f0SHong Zhang 
2084dd39f65SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_2_NaturalOrdering;
2099f5c0bcdSBarry Smith   C->ops->forwardsolve   = MatForwardSolve_SeqBAIJ_2_NaturalOrdering;
2109f5c0bcdSBarry Smith   C->ops->backwardsolve  = MatBackwardSolve_SeqBAIJ_2_NaturalOrdering;
2114dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering;
212b2b2dd24SShri Abhyankar   C->assembled           = PETSC_TRUE;
21326fbe8dcSKarl Rupp 
2149566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 2 * 2 * 2 * n)); /* from inverting diagonal blocks */
215b2b2dd24SShri Abhyankar   PetscFunctionReturn(0);
216b2b2dd24SShri Abhyankar }
217b2b2dd24SShri Abhyankar 
2189371c9d4SSatish Balay PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_inplace(Mat B, Mat A, const MatFactorInfo *info) {
219719d5645SBarry Smith   Mat             C = B;
2204eeb42bcSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
2217cdf2dbbSSatish Balay   IS              isrow = b->row, isicol = b->icol;
2225d0c19d7SBarry Smith   const PetscInt *r, *ic;
2235d0c19d7SBarry Smith   PetscInt        i, j, n = a->mbs, *bi = b->i, *bj = b->j;
224690b6cddSBarry Smith   PetscInt       *ajtmpold, *ajtmp, nz, row;
225690b6cddSBarry Smith   PetscInt       *diag_offset = b->diag, idx, *ai = a->i, *aj = a->j, *pj;
226329f5518SBarry Smith   MatScalar      *pv, *v, *rtmp, m1, m2, m3, m4, *pc, *w, *x, x1, x2, x3, x4;
2272515f8d2SSatish Balay   MatScalar       p1, p2, p3, p4;
2283f1db9ecSBarry Smith   MatScalar      *ba = b->a, *aa = a->a;
229182b8fbaSHong Zhang   PetscReal       shift = info->shiftamount;
230a455e926SHong Zhang   PetscBool       allowzeropivot, zeropivotdetected;
2314eeb42bcSBarry Smith 
2323a40ed3dSBarry Smith   PetscFunctionBegin;
2330164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
2349566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
2359566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
2369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(4 * (n + 1), &rtmp));
2374eeb42bcSBarry Smith 
2384eeb42bcSBarry Smith   for (i = 0; i < n; i++) {
2394078e994SBarry Smith     nz    = bi[i + 1] - bi[i];
2404078e994SBarry Smith     ajtmp = bj + bi[i];
2414eeb42bcSBarry Smith     for (j = 0; j < nz; j++) {
2429371c9d4SSatish Balay       x    = rtmp + 4 * ajtmp[j];
2439371c9d4SSatish Balay       x[0] = x[1] = x[2] = x[3] = 0.0;
2444eeb42bcSBarry Smith     }
2454eeb42bcSBarry Smith     /* load in initial (unfactored row) */
2464eeb42bcSBarry Smith     idx      = r[i];
2474078e994SBarry Smith     nz       = ai[idx + 1] - ai[idx];
2484078e994SBarry Smith     ajtmpold = aj + ai[idx];
2494078e994SBarry Smith     v        = aa + 4 * ai[idx];
2504eeb42bcSBarry Smith     for (j = 0; j < nz; j++) {
2514eeb42bcSBarry Smith       x    = rtmp + 4 * ic[ajtmpold[j]];
2529371c9d4SSatish Balay       x[0] = v[0];
2539371c9d4SSatish Balay       x[1] = v[1];
2549371c9d4SSatish Balay       x[2] = v[2];
2559371c9d4SSatish Balay       x[3] = v[3];
2564eeb42bcSBarry Smith       v += 4;
2574eeb42bcSBarry Smith     }
2584eeb42bcSBarry Smith     row = *ajtmp++;
2594eeb42bcSBarry Smith     while (row < i) {
2604eeb42bcSBarry Smith       pc = rtmp + 4 * row;
2619371c9d4SSatish Balay       p1 = pc[0];
2629371c9d4SSatish Balay       p2 = pc[1];
2639371c9d4SSatish Balay       p3 = pc[2];
2649371c9d4SSatish Balay       p4 = pc[3];
265d4a378daSJed Brown       if (p1 != (PetscScalar)0.0 || p2 != (PetscScalar)0.0 || p3 != (PetscScalar)0.0 || p4 != (PetscScalar)0.0) {
2664078e994SBarry Smith         pv    = ba + 4 * diag_offset[row];
2674078e994SBarry Smith         pj    = bj + diag_offset[row] + 1;
2689371c9d4SSatish Balay         x1    = pv[0];
2699371c9d4SSatish Balay         x2    = pv[1];
2709371c9d4SSatish Balay         x3    = pv[2];
2719371c9d4SSatish Balay         x4    = pv[3];
2724eeb42bcSBarry Smith         pc[0] = m1 = p1 * x1 + p3 * x2;
2734eeb42bcSBarry Smith         pc[1] = m2 = p2 * x1 + p4 * x2;
2744eeb42bcSBarry Smith         pc[2] = m3 = p1 * x3 + p3 * x4;
2754eeb42bcSBarry Smith         pc[3] = m4 = p2 * x3 + p4 * x4;
2764078e994SBarry Smith         nz         = bi[row + 1] - diag_offset[row] - 1;
2774eeb42bcSBarry Smith         pv += 4;
2784eeb42bcSBarry Smith         for (j = 0; j < nz; j++) {
2799371c9d4SSatish Balay           x1 = pv[0];
2809371c9d4SSatish Balay           x2 = pv[1];
2819371c9d4SSatish Balay           x3 = pv[2];
2829371c9d4SSatish Balay           x4 = pv[3];
2834eeb42bcSBarry Smith           x  = rtmp + 4 * pj[j];
2844eeb42bcSBarry Smith           x[0] -= m1 * x1 + m3 * x2;
2854eeb42bcSBarry Smith           x[1] -= m2 * x1 + m4 * x2;
2864eeb42bcSBarry Smith           x[2] -= m1 * x3 + m3 * x4;
2874eeb42bcSBarry Smith           x[3] -= m2 * x3 + m4 * x4;
2884eeb42bcSBarry Smith           pv += 4;
2894eeb42bcSBarry Smith         }
2909566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(16.0 * nz + 12.0));
2914eeb42bcSBarry Smith       }
2924eeb42bcSBarry Smith       row = *ajtmp++;
2934eeb42bcSBarry Smith     }
2944eeb42bcSBarry Smith     /* finished row so stick it into b->a */
2954078e994SBarry Smith     pv = ba + 4 * bi[i];
2964078e994SBarry Smith     pj = bj + bi[i];
2974078e994SBarry Smith     nz = bi[i + 1] - bi[i];
2984eeb42bcSBarry Smith     for (j = 0; j < nz; j++) {
2994eeb42bcSBarry Smith       x     = rtmp + 4 * pj[j];
3009371c9d4SSatish Balay       pv[0] = x[0];
3019371c9d4SSatish Balay       pv[1] = x[1];
3029371c9d4SSatish Balay       pv[2] = x[2];
3039371c9d4SSatish Balay       pv[3] = x[3];
3044eeb42bcSBarry Smith       pv += 4;
3054eeb42bcSBarry Smith     }
3064eeb42bcSBarry Smith     /* invert diagonal block */
3074078e994SBarry Smith     w = ba + 4 * diag_offset[i];
3089566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_2(w, shift, allowzeropivot, &zeropivotdetected));
3097b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3104eeb42bcSBarry Smith   }
3114eeb42bcSBarry Smith 
3129566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
3139566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
3149566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
31526fbe8dcSKarl Rupp 
31606e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_2_inplace;
31706e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_inplace;
3184eeb42bcSBarry Smith   C->assembled           = PETSC_TRUE;
31926fbe8dcSKarl Rupp 
3209566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 8 * b->mbs)); /* from inverting diagonal blocks */
3213a40ed3dSBarry Smith   PetscFunctionReturn(0);
3224eeb42bcSBarry Smith }
3234cd38560SBarry Smith /*
3244cd38560SBarry Smith       Version for when blocks are 2 by 2 Using natural ordering
3254cd38560SBarry Smith */
3269371c9d4SSatish Balay PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace(Mat C, Mat A, const MatFactorInfo *info) {
3274cd38560SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
328690b6cddSBarry Smith   PetscInt     i, j, n = a->mbs, *bi = b->i, *bj = b->j;
329690b6cddSBarry Smith   PetscInt    *ajtmpold, *ajtmp, nz, row;
330690b6cddSBarry Smith   PetscInt    *diag_offset = b->diag, *ai = a->i, *aj = a->j, *pj;
331329f5518SBarry Smith   MatScalar   *pv, *v, *rtmp, *pc, *w, *x;
3322515f8d2SSatish Balay   MatScalar    p1, p2, p3, p4, m1, m2, m3, m4, x1, x2, x3, x4;
3334cd38560SBarry Smith   MatScalar   *ba = b->a, *aa = a->a;
334182b8fbaSHong Zhang   PetscReal    shift = info->shiftamount;
335a455e926SHong Zhang   PetscBool    allowzeropivot, zeropivotdetected;
3364cd38560SBarry Smith 
3374cd38560SBarry Smith   PetscFunctionBegin;
3380164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
3399566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(4 * (n + 1), &rtmp));
3404cd38560SBarry Smith   for (i = 0; i < n; i++) {
3414cd38560SBarry Smith     nz    = bi[i + 1] - bi[i];
3424cd38560SBarry Smith     ajtmp = bj + bi[i];
3434cd38560SBarry Smith     for (j = 0; j < nz; j++) {
3444cd38560SBarry Smith       x    = rtmp + 4 * ajtmp[j];
3454cd38560SBarry Smith       x[0] = x[1] = x[2] = x[3] = 0.0;
3464cd38560SBarry Smith     }
3474cd38560SBarry Smith     /* load in initial (unfactored row) */
3484cd38560SBarry Smith     nz       = ai[i + 1] - ai[i];
3494cd38560SBarry Smith     ajtmpold = aj + ai[i];
3504cd38560SBarry Smith     v        = aa + 4 * ai[i];
3514cd38560SBarry Smith     for (j = 0; j < nz; j++) {
3524cd38560SBarry Smith       x    = rtmp + 4 * ajtmpold[j];
3539371c9d4SSatish Balay       x[0] = v[0];
3549371c9d4SSatish Balay       x[1] = v[1];
3559371c9d4SSatish Balay       x[2] = v[2];
3569371c9d4SSatish Balay       x[3] = v[3];
3574cd38560SBarry Smith       v += 4;
3584cd38560SBarry Smith     }
3594cd38560SBarry Smith     row = *ajtmp++;
3604cd38560SBarry Smith     while (row < i) {
3614cd38560SBarry Smith       pc = rtmp + 4 * row;
3629371c9d4SSatish Balay       p1 = pc[0];
3639371c9d4SSatish Balay       p2 = pc[1];
3649371c9d4SSatish Balay       p3 = pc[2];
3659371c9d4SSatish Balay       p4 = pc[3];
366d4a378daSJed Brown       if (p1 != (PetscScalar)0.0 || p2 != (PetscScalar)0.0 || p3 != (PetscScalar)0.0 || p4 != (PetscScalar)0.0) {
3674cd38560SBarry Smith         pv    = ba + 4 * diag_offset[row];
3684cd38560SBarry Smith         pj    = bj + diag_offset[row] + 1;
3699371c9d4SSatish Balay         x1    = pv[0];
3709371c9d4SSatish Balay         x2    = pv[1];
3719371c9d4SSatish Balay         x3    = pv[2];
3729371c9d4SSatish Balay         x4    = pv[3];
3734cd38560SBarry Smith         pc[0] = m1 = p1 * x1 + p3 * x2;
3744cd38560SBarry Smith         pc[1] = m2 = p2 * x1 + p4 * x2;
3754cd38560SBarry Smith         pc[2] = m3 = p1 * x3 + p3 * x4;
3764cd38560SBarry Smith         pc[3] = m4 = p2 * x3 + p4 * x4;
3774cd38560SBarry Smith         nz         = bi[row + 1] - diag_offset[row] - 1;
3784cd38560SBarry Smith         pv += 4;
3794cd38560SBarry Smith         for (j = 0; j < nz; j++) {
3809371c9d4SSatish Balay           x1 = pv[0];
3819371c9d4SSatish Balay           x2 = pv[1];
3829371c9d4SSatish Balay           x3 = pv[2];
3839371c9d4SSatish Balay           x4 = pv[3];
3844cd38560SBarry Smith           x  = rtmp + 4 * pj[j];
3854cd38560SBarry Smith           x[0] -= m1 * x1 + m3 * x2;
3864cd38560SBarry Smith           x[1] -= m2 * x1 + m4 * x2;
3874cd38560SBarry Smith           x[2] -= m1 * x3 + m3 * x4;
3884cd38560SBarry Smith           x[3] -= m2 * x3 + m4 * x4;
3894cd38560SBarry Smith           pv += 4;
3904cd38560SBarry Smith         }
3919566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(16.0 * nz + 12.0));
3924cd38560SBarry Smith       }
3934cd38560SBarry Smith       row = *ajtmp++;
3944cd38560SBarry Smith     }
3954cd38560SBarry Smith     /* finished row so stick it into b->a */
3964cd38560SBarry Smith     pv = ba + 4 * bi[i];
3974cd38560SBarry Smith     pj = bj + bi[i];
3984cd38560SBarry Smith     nz = bi[i + 1] - bi[i];
3994cd38560SBarry Smith     for (j = 0; j < nz; j++) {
4004cd38560SBarry Smith       x     = rtmp + 4 * pj[j];
4019371c9d4SSatish Balay       pv[0] = x[0];
4029371c9d4SSatish Balay       pv[1] = x[1];
4039371c9d4SSatish Balay       pv[2] = x[2];
4049371c9d4SSatish Balay       pv[3] = x[3];
4052efa7f71SHong Zhang       /*
4062efa7f71SHong Zhang       printf(" col %d:",pj[j]);
4072efa7f71SHong Zhang       PetscInt j1;
4082efa7f71SHong Zhang       for (j1=0; j1<4; j1++) printf(" %g,",*(pv+j1));
4092efa7f71SHong Zhang       printf("\n");
4102efa7f71SHong Zhang       */
4114cd38560SBarry Smith       pv += 4;
4124cd38560SBarry Smith     }
4134cd38560SBarry Smith     /* invert diagonal block */
4144cd38560SBarry Smith     w = ba + 4 * diag_offset[i];
4159566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_2(w, shift, allowzeropivot, &zeropivotdetected));
4167b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
4174cd38560SBarry Smith   }
4184cd38560SBarry Smith 
4199566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
42026fbe8dcSKarl Rupp 
42106e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_2_NaturalOrdering_inplace;
42206e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering_inplace;
4234cd38560SBarry Smith   C->assembled           = PETSC_TRUE;
42426fbe8dcSKarl Rupp 
4259566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 8 * b->mbs)); /* from inverting diagonal blocks */
4264cd38560SBarry Smith   PetscFunctionReturn(0);
4274cd38560SBarry Smith }
4287fc0212eSBarry Smith 
4297fc0212eSBarry Smith /* ----------------------------------------------------------- */
4307fc0212eSBarry Smith /*
4317fc0212eSBarry Smith      Version for when blocks are 1 by 1.
4327fc0212eSBarry Smith */
4339371c9d4SSatish Balay PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat B, Mat A, const MatFactorInfo *info) {
434048b5e81SShri Abhyankar   Mat              C = B;
435048b5e81SShri Abhyankar   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
436048b5e81SShri Abhyankar   IS               isrow = b->row, isicol = b->icol;
437048b5e81SShri Abhyankar   const PetscInt  *r, *ic, *ics;
438048b5e81SShri Abhyankar   const PetscInt   n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j, *bdiag = b->diag;
439048b5e81SShri Abhyankar   PetscInt         i, j, k, nz, nzL, row, *pj;
440048b5e81SShri Abhyankar   const PetscInt  *ajtmp, *bjtmp;
441048b5e81SShri Abhyankar   MatScalar       *rtmp, *pc, multiplier, *pv;
442048b5e81SShri Abhyankar   const MatScalar *aa = a->a, *v;
443ace3abfcSBarry Smith   PetscBool        row_identity, col_identity;
444048b5e81SShri Abhyankar   FactorShiftCtx   sctx;
445048b5e81SShri Abhyankar   const PetscInt  *ddiag;
446048b5e81SShri Abhyankar   PetscReal        rs;
447048b5e81SShri Abhyankar   MatScalar        d;
448048b5e81SShri Abhyankar 
449048b5e81SShri Abhyankar   PetscFunctionBegin;
450048b5e81SShri Abhyankar   /* MatPivotSetUp(): initialize shift context sctx */
4519566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&sctx, sizeof(FactorShiftCtx)));
452048b5e81SShri Abhyankar 
453048b5e81SShri Abhyankar   if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
454048b5e81SShri Abhyankar     ddiag          = a->diag;
455048b5e81SShri Abhyankar     sctx.shift_top = info->zeropivot;
456048b5e81SShri Abhyankar     for (i = 0; i < n; i++) {
457048b5e81SShri Abhyankar       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
458048b5e81SShri Abhyankar       d  = (aa)[ddiag[i]];
459048b5e81SShri Abhyankar       rs = -PetscAbsScalar(d) - PetscRealPart(d);
460048b5e81SShri Abhyankar       v  = aa + ai[i];
461048b5e81SShri Abhyankar       nz = ai[i + 1] - ai[i];
46226fbe8dcSKarl Rupp       for (j = 0; j < nz; j++) rs += PetscAbsScalar(v[j]);
463048b5e81SShri Abhyankar       if (rs > sctx.shift_top) sctx.shift_top = rs;
464048b5e81SShri Abhyankar     }
465048b5e81SShri Abhyankar     sctx.shift_top *= 1.1;
466048b5e81SShri Abhyankar     sctx.nshift_max = 5;
467048b5e81SShri Abhyankar     sctx.shift_lo   = 0.;
468048b5e81SShri Abhyankar     sctx.shift_hi   = 1.;
469048b5e81SShri Abhyankar   }
470048b5e81SShri Abhyankar 
4719566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
4729566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
4739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n + 1, &rtmp));
474048b5e81SShri Abhyankar   ics = ic;
475048b5e81SShri Abhyankar 
476048b5e81SShri Abhyankar   do {
477048b5e81SShri Abhyankar     sctx.newshift = PETSC_FALSE;
478048b5e81SShri Abhyankar     for (i = 0; i < n; i++) {
479048b5e81SShri Abhyankar       /* zero rtmp */
480048b5e81SShri Abhyankar       /* L part */
481048b5e81SShri Abhyankar       nz    = bi[i + 1] - bi[i];
482048b5e81SShri Abhyankar       bjtmp = bj + bi[i];
483048b5e81SShri Abhyankar       for (j = 0; j < nz; j++) rtmp[bjtmp[j]] = 0.0;
484048b5e81SShri Abhyankar 
485048b5e81SShri Abhyankar       /* U part */
486048b5e81SShri Abhyankar       nz    = bdiag[i] - bdiag[i + 1];
487048b5e81SShri Abhyankar       bjtmp = bj + bdiag[i + 1] + 1;
488048b5e81SShri Abhyankar       for (j = 0; j < nz; j++) rtmp[bjtmp[j]] = 0.0;
489048b5e81SShri Abhyankar 
490048b5e81SShri Abhyankar       /* load in initial (unfactored row) */
491048b5e81SShri Abhyankar       nz    = ai[r[i] + 1] - ai[r[i]];
492048b5e81SShri Abhyankar       ajtmp = aj + ai[r[i]];
493048b5e81SShri Abhyankar       v     = aa + ai[r[i]];
49426fbe8dcSKarl Rupp       for (j = 0; j < nz; j++) rtmp[ics[ajtmp[j]]] = v[j];
49526fbe8dcSKarl Rupp 
496048b5e81SShri Abhyankar       /* ZeropivotApply() */
497048b5e81SShri Abhyankar       rtmp[i] += sctx.shift_amount; /* shift the diagonal of the matrix */
498048b5e81SShri Abhyankar 
499048b5e81SShri Abhyankar       /* elimination */
500048b5e81SShri Abhyankar       bjtmp = bj + bi[i];
501048b5e81SShri Abhyankar       row   = *bjtmp++;
502048b5e81SShri Abhyankar       nzL   = bi[i + 1] - bi[i];
503048b5e81SShri Abhyankar       for (k = 0; k < nzL; k++) {
504048b5e81SShri Abhyankar         pc = rtmp + row;
505d4a378daSJed Brown         if (*pc != (PetscScalar)0.0) {
506048b5e81SShri Abhyankar           pv         = b->a + bdiag[row];
507048b5e81SShri Abhyankar           multiplier = *pc * (*pv);
508048b5e81SShri Abhyankar           *pc        = multiplier;
50926fbe8dcSKarl Rupp 
510048b5e81SShri Abhyankar           pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */
511048b5e81SShri Abhyankar           pv = b->a + bdiag[row + 1] + 1;
512048b5e81SShri Abhyankar           nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries in U(row,:) excluding diag */
513048b5e81SShri Abhyankar           for (j = 0; j < nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
5149566063dSJacob Faibussowitsch           PetscCall(PetscLogFlops(2.0 * nz));
515048b5e81SShri Abhyankar         }
516048b5e81SShri Abhyankar         row = *bjtmp++;
517048b5e81SShri Abhyankar       }
518048b5e81SShri Abhyankar 
519048b5e81SShri Abhyankar       /* finished row so stick it into b->a */
520048b5e81SShri Abhyankar       rs = 0.0;
521048b5e81SShri Abhyankar       /* L part */
522048b5e81SShri Abhyankar       pv = b->a + bi[i];
523048b5e81SShri Abhyankar       pj = b->j + bi[i];
524048b5e81SShri Abhyankar       nz = bi[i + 1] - bi[i];
525048b5e81SShri Abhyankar       for (j = 0; j < nz; j++) {
5269371c9d4SSatish Balay         pv[j] = rtmp[pj[j]];
5279371c9d4SSatish Balay         rs += PetscAbsScalar(pv[j]);
528048b5e81SShri Abhyankar       }
529048b5e81SShri Abhyankar 
530048b5e81SShri Abhyankar       /* U part */
531048b5e81SShri Abhyankar       pv = b->a + bdiag[i + 1] + 1;
532048b5e81SShri Abhyankar       pj = b->j + bdiag[i + 1] + 1;
533048b5e81SShri Abhyankar       nz = bdiag[i] - bdiag[i + 1] - 1;
534048b5e81SShri Abhyankar       for (j = 0; j < nz; j++) {
5359371c9d4SSatish Balay         pv[j] = rtmp[pj[j]];
5369371c9d4SSatish Balay         rs += PetscAbsScalar(pv[j]);
537048b5e81SShri Abhyankar       }
538048b5e81SShri Abhyankar 
539048b5e81SShri Abhyankar       sctx.rs = rs;
540048b5e81SShri Abhyankar       sctx.pv = rtmp[i];
5419566063dSJacob Faibussowitsch       PetscCall(MatPivotCheck(B, A, info, &sctx, i));
542048b5e81SShri Abhyankar       if (sctx.newshift) break; /* break for-loop */
543048b5e81SShri Abhyankar       rtmp[i] = sctx.pv;        /* sctx.pv might be updated in the case of MAT_SHIFT_INBLOCKS */
544048b5e81SShri Abhyankar 
545a5b23f4aSJose E. Roman       /* Mark diagonal and invert diagonal for simpler triangular solves */
546048b5e81SShri Abhyankar       pv  = b->a + bdiag[i];
547d4a378daSJed Brown       *pv = (PetscScalar)1.0 / rtmp[i];
548048b5e81SShri Abhyankar 
549048b5e81SShri Abhyankar     } /* endof for (i=0; i<n; i++) { */
550048b5e81SShri Abhyankar 
551048b5e81SShri Abhyankar     /* MatPivotRefine() */
552048b5e81SShri Abhyankar     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE && !sctx.newshift && sctx.shift_fraction > 0 && sctx.nshift < sctx.nshift_max) {
553048b5e81SShri Abhyankar       /*
554048b5e81SShri Abhyankar        * if no shift in this attempt & shifting & started shifting & can refine,
555048b5e81SShri Abhyankar        * then try lower shift
556048b5e81SShri Abhyankar        */
557048b5e81SShri Abhyankar       sctx.shift_hi       = sctx.shift_fraction;
558048b5e81SShri Abhyankar       sctx.shift_fraction = (sctx.shift_hi + sctx.shift_lo) / 2.;
559048b5e81SShri Abhyankar       sctx.shift_amount   = sctx.shift_fraction * sctx.shift_top;
560048b5e81SShri Abhyankar       sctx.newshift       = PETSC_TRUE;
561048b5e81SShri Abhyankar       sctx.nshift++;
562048b5e81SShri Abhyankar     }
563048b5e81SShri Abhyankar   } while (sctx.newshift);
564048b5e81SShri Abhyankar 
5659566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
5669566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
5679566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
568048b5e81SShri Abhyankar 
5699566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isrow, &row_identity));
5709566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isicol, &col_identity));
571048b5e81SShri Abhyankar   if (row_identity && col_identity) {
572048b5e81SShri Abhyankar     C->ops->solve          = MatSolve_SeqBAIJ_1_NaturalOrdering;
5739f5c0bcdSBarry Smith     C->ops->forwardsolve   = MatForwardSolve_SeqBAIJ_1_NaturalOrdering;
5749f5c0bcdSBarry Smith     C->ops->backwardsolve  = MatBackwardSolve_SeqBAIJ_1_NaturalOrdering;
57593fd935bSShri Abhyankar     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering;
576048b5e81SShri Abhyankar   } else {
577048b5e81SShri Abhyankar     C->ops->solve          = MatSolve_SeqBAIJ_1;
57893fd935bSShri Abhyankar     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1;
579048b5e81SShri Abhyankar   }
580048b5e81SShri Abhyankar   C->assembled = PETSC_TRUE;
5819566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(C->cmap->n));
582048b5e81SShri Abhyankar 
583048b5e81SShri Abhyankar   /* MatShiftView(A,info,&sctx) */
584048b5e81SShri Abhyankar   if (sctx.nshift) {
585048b5e81SShri Abhyankar     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
5869566063dSJacob 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));
587048b5e81SShri Abhyankar     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
5889566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shift_nz tries %" PetscInt_FMT ", shift_amount %g\n", sctx.nshift, (double)sctx.shift_amount));
589048b5e81SShri Abhyankar     } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) {
5909566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shift_inblocks applied %" PetscInt_FMT ", each shift_amount %g\n", sctx.nshift, (double)info->shiftamount));
591048b5e81SShri Abhyankar     }
592048b5e81SShri Abhyankar   }
593048b5e81SShri Abhyankar   PetscFunctionReturn(0);
594048b5e81SShri Abhyankar }
595048b5e81SShri Abhyankar 
5969371c9d4SSatish Balay PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1_inplace(Mat C, Mat A, const MatFactorInfo *info) {
5977fc0212eSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
5987cdf2dbbSSatish Balay   IS              isrow = b->row, isicol = b->icol;
5995d0c19d7SBarry Smith   const PetscInt *r, *ic;
6005d0c19d7SBarry Smith   PetscInt        i, j, n = a->mbs, *bi = b->i, *bj = b->j;
601690b6cddSBarry Smith   PetscInt       *ajtmpold, *ajtmp, nz, row, *ai = a->i, *aj = a->j;
602690b6cddSBarry Smith   PetscInt       *diag_offset = b->diag, diag, *pj;
603329f5518SBarry Smith   MatScalar      *pv, *v, *rtmp, multiplier, *pc;
6043f1db9ecSBarry Smith   MatScalar      *ba = b->a, *aa = a->a;
605ace3abfcSBarry Smith   PetscBool       row_identity, col_identity;
6067fc0212eSBarry Smith 
6073a40ed3dSBarry Smith   PetscFunctionBegin;
6089566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
6099566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
6109566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n + 1, &rtmp));
6117fc0212eSBarry Smith 
6127fc0212eSBarry Smith   for (i = 0; i < n; i++) {
6134078e994SBarry Smith     nz    = bi[i + 1] - bi[i];
6144078e994SBarry Smith     ajtmp = bj + bi[i];
6157fc0212eSBarry Smith     for (j = 0; j < nz; j++) rtmp[ajtmp[j]] = 0.0;
6167fc0212eSBarry Smith 
6177fc0212eSBarry Smith     /* load in initial (unfactored row) */
6184078e994SBarry Smith     nz       = ai[r[i] + 1] - ai[r[i]];
6194078e994SBarry Smith     ajtmpold = aj + ai[r[i]];
6204078e994SBarry Smith     v        = aa + ai[r[i]];
6217fc0212eSBarry Smith     for (j = 0; j < nz; j++) rtmp[ic[ajtmpold[j]]] = v[j];
6227fc0212eSBarry Smith 
6237fc0212eSBarry Smith     row = *ajtmp++;
6247fc0212eSBarry Smith     while (row < i) {
6257fc0212eSBarry Smith       pc = rtmp + row;
6267fc0212eSBarry Smith       if (*pc != 0.0) {
6274078e994SBarry Smith         pv         = ba + diag_offset[row];
6284078e994SBarry Smith         pj         = bj + diag_offset[row] + 1;
6297fc0212eSBarry Smith         multiplier = *pc * *pv++;
6307fc0212eSBarry Smith         *pc        = multiplier;
6314078e994SBarry Smith         nz         = bi[row + 1] - diag_offset[row] - 1;
6327fc0212eSBarry Smith         for (j = 0; j < nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
6339566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(1.0 + 2.0 * nz));
6347fc0212eSBarry Smith       }
6357fc0212eSBarry Smith       row = *ajtmp++;
6367fc0212eSBarry Smith     }
6377fc0212eSBarry Smith     /* finished row so stick it into b->a */
6384078e994SBarry Smith     pv = ba + bi[i];
6394078e994SBarry Smith     pj = bj + bi[i];
6404078e994SBarry Smith     nz = bi[i + 1] - bi[i];
64126fbe8dcSKarl Rupp     for (j = 0; j < nz; j++) pv[j] = rtmp[pj[j]];
6424078e994SBarry Smith     diag = diag_offset[i] - bi[i];
6437fc0212eSBarry Smith     /* check pivot entry for current row */
64408401ef6SPierre Jolivet     PetscCheck(pv[diag] != 0.0, PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot: row in original ordering %" PetscInt_FMT " in permuted ordering %" PetscInt_FMT, r[i], i);
6457fc0212eSBarry Smith     pv[diag] = 1.0 / pv[diag];
6467fc0212eSBarry Smith   }
6477fc0212eSBarry Smith 
6489566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
6499566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
6509566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
6519566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isrow, &row_identity));
6529566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isicol, &col_identity));
653f58c8c74SBarry Smith   if (row_identity && col_identity) {
65406e38f1dSHong Zhang     C->ops->solve          = MatSolve_SeqBAIJ_1_NaturalOrdering_inplace;
65506e38f1dSHong Zhang     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering_inplace;
656f58c8c74SBarry Smith   } else {
65706e38f1dSHong Zhang     C->ops->solve          = MatSolve_SeqBAIJ_1_inplace;
65806e38f1dSHong Zhang     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_inplace;
659f58c8c74SBarry Smith   }
6607fc0212eSBarry Smith   C->assembled = PETSC_TRUE;
6619566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(C->cmap->n));
6623a40ed3dSBarry Smith   PetscFunctionReturn(0);
6637fc0212eSBarry Smith }
6647fc0212eSBarry Smith 
6659371c9d4SSatish Balay static PetscErrorCode MatFactorGetSolverType_petsc(Mat A, MatSolverType *type) {
6664ac6704cSBarry Smith   PetscFunctionBegin;
6674ac6704cSBarry Smith   *type = MATSOLVERPETSC;
6684ac6704cSBarry Smith   PetscFunctionReturn(0);
6694ac6704cSBarry Smith }
6704ac6704cSBarry Smith 
6719371c9d4SSatish Balay PETSC_INTERN PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A, MatFactorType ftype, Mat *B) {
672d0f46423SBarry Smith   PetscInt n = A->rmap->n;
673b24902e0SBarry Smith 
674b24902e0SBarry Smith   PetscFunctionBegin;
675b499a2aeSHong Zhang #if defined(PETSC_USE_COMPLEX)
676b94d7dedSBarry Smith   PetscCheck(A->hermitian != PETSC_BOOL3_TRUE || !(ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC), PETSC_COMM_SELF, PETSC_ERR_SUP, "Hermitian Factor is not supported");
677b499a2aeSHong Zhang #endif
6789566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
6799566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*B, n, n, n, n));
680c8342467SHong Zhang   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) {
6819566063dSJacob Faibussowitsch     PetscCall(MatSetType(*B, MATSEQBAIJ));
68226fbe8dcSKarl Rupp 
6834dd39f65SShri Abhyankar     (*B)->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqBAIJ;
6844dd39f65SShri Abhyankar     (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqBAIJ;
6859566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGND, (char **)&(*B)->preferredordering[MAT_FACTOR_LU]));
6869566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGNATURAL, (char **)&(*B)->preferredordering[MAT_FACTOR_ILU]));
6879566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGNATURAL, (char **)&(*B)->preferredordering[MAT_FACTOR_ILUDT]));
688b24902e0SBarry Smith   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
6899566063dSJacob Faibussowitsch     PetscCall(MatSetType(*B, MATSEQSBAIJ));
6909566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(*B, A->rmap->bs, MAT_SKIP_ALLOCATION, NULL));
69126fbe8dcSKarl Rupp 
6925c9eb25fSBarry Smith     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqBAIJ;
6935c9eb25fSBarry Smith     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ;
6944ac6704cSBarry Smith     /*  Future optimization would be direct symbolic and numerical factorization for BAIJ to support orderings and Cholesky, instead of first converting to SBAIJ */
6959566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGNATURAL, (char **)&(*B)->preferredordering[MAT_FACTOR_CHOLESKY]));
6969566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGNATURAL, (char **)&(*B)->preferredordering[MAT_FACTOR_ICC]));
697e32f2f54SBarry Smith   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported");
698d5f3da31SBarry Smith   (*B)->factortype     = ftype;
699f73b0415SBarry Smith   (*B)->canuseordering = PETSC_TRUE;
70000c67f3bSHong Zhang 
7019566063dSJacob Faibussowitsch   PetscCall(PetscFree((*B)->solvertype));
7029566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERPETSC, &(*B)->solvertype));
7039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)*B, "MatFactorGetSolverType_C", MatFactorGetSolverType_petsc));
704b24902e0SBarry Smith   PetscFunctionReturn(0);
705b24902e0SBarry Smith }
706273d9f13SBarry Smith 
7075d517e7eSBarry Smith /* ----------------------------------------------------------- */
7089371c9d4SSatish Balay PetscErrorCode MatLUFactor_SeqBAIJ(Mat A, IS row, IS col, const MatFactorInfo *info) {
7095d517e7eSBarry Smith   Mat C;
7105d517e7eSBarry Smith 
7113a40ed3dSBarry Smith   PetscFunctionBegin;
7129566063dSJacob Faibussowitsch   PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_LU, &C));
7139566063dSJacob Faibussowitsch   PetscCall(MatLUFactorSymbolic(C, A, row, col, info));
7149566063dSJacob Faibussowitsch   PetscCall(MatLUFactorNumeric(C, A, info));
71526fbe8dcSKarl Rupp 
716db4efbfdSBarry Smith   A->ops->solve          = C->ops->solve;
717db4efbfdSBarry Smith   A->ops->solvetranspose = C->ops->solvetranspose;
71826fbe8dcSKarl Rupp 
7199566063dSJacob Faibussowitsch   PetscCall(MatHeaderMerge(A, &C));
7209566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)A, (PetscObject)((Mat_SeqBAIJ *)(A->data))->icol));
7213a40ed3dSBarry Smith   PetscFunctionReturn(0);
7225d517e7eSBarry Smith }
7234cd38560SBarry Smith 
724c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h>
7259371c9d4SSatish Balay PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat C, Mat A, const MatFactorInfo *info) {
72678910aadSHong Zhang   Mat_SeqBAIJ    *a  = (Mat_SeqBAIJ *)A->data;
72778910aadSHong Zhang   Mat_SeqSBAIJ   *b  = (Mat_SeqSBAIJ *)C->data;
72878910aadSHong Zhang   IS              ip = b->row;
7295d0c19d7SBarry Smith   const PetscInt *rip;
7305d0c19d7SBarry Smith   PetscInt        i, j, mbs = a->mbs, bs = A->rmap->bs, *bi = b->i, *bj = b->j, *bcol;
73178910aadSHong Zhang   PetscInt       *ai = a->i, *aj = a->j;
73278910aadSHong Zhang   PetscInt        k, jmin, jmax, *jl, *il, col, nexti, ili, nz;
73378910aadSHong Zhang   MatScalar      *rtmp, *ba = b->a, *bval, *aa = a->a, dk, uikdi;
73407b50cabSHong Zhang   PetscReal       rs;
7350e95ead3SHong Zhang   FactorShiftCtx  sctx;
73678910aadSHong Zhang 
737c05c3958SHong Zhang   PetscFunctionBegin;
73807b50cabSHong Zhang   if (bs > 1) { /* convert A to a SBAIJ matrix and apply Cholesky factorization from it */
739*48a46eb9SPierre Jolivet     if (!a->sbaijMat) PetscCall(MatConvert(A, MATSEQSBAIJ, MAT_INITIAL_MATRIX, &a->sbaijMat));
7409566063dSJacob Faibussowitsch     PetscCall((a->sbaijMat)->ops->choleskyfactornumeric(C, a->sbaijMat, info));
7419566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&a->sbaijMat));
7426ad2eaddSHong Zhang     PetscFunctionReturn(0);
7436ad2eaddSHong Zhang   }
74478910aadSHong Zhang 
74507b50cabSHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
7469566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&sctx, sizeof(FactorShiftCtx)));
7473cea4cbeSHong Zhang 
7489566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(ip, &rip));
7499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(mbs, &rtmp, mbs, &il, mbs, &jl));
75078910aadSHong Zhang 
75175567043SBarry Smith   sctx.shift_amount = 0.;
7523cea4cbeSHong Zhang   sctx.nshift       = 0;
75378910aadSHong Zhang   do {
75407b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
75578910aadSHong Zhang     for (i = 0; i < mbs; i++) {
7569371c9d4SSatish Balay       rtmp[i] = 0.0;
7579371c9d4SSatish Balay       jl[i]   = mbs;
7589371c9d4SSatish Balay       il[0]   = 0;
75978910aadSHong Zhang     }
76078910aadSHong Zhang 
76178910aadSHong Zhang     for (k = 0; k < mbs; k++) {
76278910aadSHong Zhang       bval = ba + bi[k];
76378910aadSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
7649371c9d4SSatish Balay       jmin = ai[rip[k]];
7659371c9d4SSatish Balay       jmax = ai[rip[k] + 1];
76678910aadSHong Zhang       for (j = jmin; j < jmax; j++) {
76778910aadSHong Zhang         col = rip[aj[j]];
76878910aadSHong Zhang         if (col >= k) { /* only take upper triangular entry */
76978910aadSHong Zhang           rtmp[col] = aa[j];
77078910aadSHong Zhang           *bval++   = 0.0; /* for in-place factorization */
77178910aadSHong Zhang         }
77278910aadSHong Zhang       }
7733cea4cbeSHong Zhang 
7743cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
7753cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
77678910aadSHong Zhang 
77778910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
77878910aadSHong Zhang       dk = rtmp[k];
77978910aadSHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
78078910aadSHong Zhang 
78178910aadSHong Zhang       while (i < k) {
78278910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
78378910aadSHong Zhang 
78478910aadSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
78578910aadSHong Zhang         ili   = il[i];                /* index of first nonzero element in U(i,k:bms-1) */
78678910aadSHong Zhang         uikdi = -ba[ili] * ba[bi[i]]; /* diagonal(k) */
78778910aadSHong Zhang         dk += uikdi * ba[ili];
78878910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
78978910aadSHong Zhang 
79078910aadSHong Zhang         /* add multiple of row i to k-th row */
7919371c9d4SSatish Balay         jmin = ili + 1;
7929371c9d4SSatish Balay         jmax = bi[i + 1];
79378910aadSHong Zhang         if (jmin < jmax) {
79478910aadSHong Zhang           for (j = jmin; j < jmax; j++) rtmp[bj[j]] += uikdi * ba[j];
79578910aadSHong Zhang           /* update il and jl for row i */
79678910aadSHong Zhang           il[i] = jmin;
7979371c9d4SSatish Balay           j     = bj[jmin];
7989371c9d4SSatish Balay           jl[i] = jl[j];
7999371c9d4SSatish Balay           jl[j] = i;
80078910aadSHong Zhang         }
80178910aadSHong Zhang         i = nexti;
80278910aadSHong Zhang       }
80378910aadSHong Zhang 
8043cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
8053cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
8063cea4cbeSHong Zhang       rs   = 0.0;
80778910aadSHong Zhang       jmin = bi[k] + 1;
80878910aadSHong Zhang       nz   = bi[k + 1] - jmin;
80978910aadSHong Zhang       if (nz) {
81078910aadSHong Zhang         bcol = bj + jmin;
81178910aadSHong Zhang         while (nz--) {
81289f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
81389f845c8SHong Zhang           bcol++;
81478910aadSHong Zhang         }
81578910aadSHong Zhang       }
8163cea4cbeSHong Zhang 
8173cea4cbeSHong Zhang       sctx.rs = rs;
8183cea4cbeSHong Zhang       sctx.pv = dk;
8199566063dSJacob Faibussowitsch       PetscCall(MatPivotCheck(C, A, info, &sctx, k));
82007b50cabSHong Zhang       if (sctx.newshift) break;
8210e95ead3SHong Zhang       dk = sctx.pv;
82278910aadSHong Zhang 
82378910aadSHong Zhang       /* copy data into U(k,:) */
82478910aadSHong Zhang       ba[bi[k]] = 1.0 / dk; /* U(k,k) */
8259371c9d4SSatish Balay       jmin      = bi[k] + 1;
8269371c9d4SSatish Balay       jmax      = bi[k + 1];
82778910aadSHong Zhang       if (jmin < jmax) {
82878910aadSHong Zhang         for (j = jmin; j < jmax; j++) {
8299371c9d4SSatish Balay           col       = bj[j];
8309371c9d4SSatish Balay           ba[j]     = rtmp[col];
8319371c9d4SSatish Balay           rtmp[col] = 0.0;
83278910aadSHong Zhang         }
83378910aadSHong Zhang         /* add the k-th row into il and jl */
83478910aadSHong Zhang         il[k] = jmin;
8359371c9d4SSatish Balay         i     = bj[jmin];
8369371c9d4SSatish Balay         jl[k] = jl[i];
8379371c9d4SSatish Balay         jl[i] = k;
83878910aadSHong Zhang       }
83978910aadSHong Zhang     }
84007b50cabSHong Zhang   } while (sctx.newshift);
8419566063dSJacob Faibussowitsch   PetscCall(PetscFree3(rtmp, il, jl));
84278910aadSHong Zhang 
8439566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(ip, &rip));
84426fbe8dcSKarl Rupp 
84578910aadSHong Zhang   C->assembled    = PETSC_TRUE;
84678910aadSHong Zhang   C->preallocated = PETSC_TRUE;
84726fbe8dcSKarl Rupp 
8489566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(C->rmap->N));
8493cea4cbeSHong Zhang   if (sctx.nshift) {
850f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
8519566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shiftpd tries %" PetscInt_FMT ", shift_amount %g\n", sctx.nshift, (double)sctx.shift_amount));
852f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
8539566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shiftnz tries %" PetscInt_FMT ", shift_amount %g\n", sctx.nshift, (double)sctx.shift_amount));
85478910aadSHong Zhang     }
85578910aadSHong Zhang   }
856c05c3958SHong Zhang   PetscFunctionReturn(0);
857c05c3958SHong Zhang }
8584cd38560SBarry Smith 
8599371c9d4SSatish Balay PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat C, Mat A, const MatFactorInfo *info) {
86078910aadSHong Zhang   Mat_SeqBAIJ   *a = (Mat_SeqBAIJ *)A->data;
86178910aadSHong Zhang   Mat_SeqSBAIJ  *b = (Mat_SeqSBAIJ *)C->data;
8623cea4cbeSHong Zhang   PetscInt       i, j, am = a->mbs;
86378910aadSHong Zhang   PetscInt      *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j;
8643cea4cbeSHong Zhang   PetscInt       k, jmin, *jl, *il, nexti, ili, *acol, *bcol, nz;
86578910aadSHong Zhang   MatScalar     *rtmp, *ba = b->a, *aa = a->a, dk, uikdi, *aval, *bval;
8660e95ead3SHong Zhang   PetscReal      rs;
8670e95ead3SHong Zhang   FactorShiftCtx sctx;
86878910aadSHong Zhang 
869c05c3958SHong Zhang   PetscFunctionBegin;
8700e95ead3SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
8719566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&sctx, sizeof(FactorShiftCtx)));
8723cea4cbeSHong Zhang 
8739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(am, &rtmp, am, &il, am, &jl));
87478910aadSHong Zhang 
87578910aadSHong Zhang   do {
87607b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
87778910aadSHong Zhang     for (i = 0; i < am; i++) {
8789371c9d4SSatish Balay       rtmp[i] = 0.0;
8799371c9d4SSatish Balay       jl[i]   = am;
8809371c9d4SSatish Balay       il[0]   = 0;
88178910aadSHong Zhang     }
88278910aadSHong Zhang 
88378910aadSHong Zhang     for (k = 0; k < am; k++) {
88478910aadSHong Zhang       /* initialize k-th row with elements nonzero in row perm(k) of A */
88578910aadSHong Zhang       nz   = ai[k + 1] - ai[k];
88678910aadSHong Zhang       acol = aj + ai[k];
88778910aadSHong Zhang       aval = aa + ai[k];
88878910aadSHong Zhang       bval = ba + bi[k];
88978910aadSHong Zhang       while (nz--) {
89078910aadSHong Zhang         if (*acol < k) { /* skip lower triangular entries */
8919371c9d4SSatish Balay           acol++;
8929371c9d4SSatish Balay           aval++;
89378910aadSHong Zhang         } else {
89478910aadSHong Zhang           rtmp[*acol++] = *aval++;
89578910aadSHong Zhang           *bval++       = 0.0; /* for in-place factorization */
89678910aadSHong Zhang         }
89778910aadSHong Zhang       }
8983cea4cbeSHong Zhang 
8993cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
9003cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
90178910aadSHong Zhang 
90278910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
90378910aadSHong Zhang       dk = rtmp[k];
90478910aadSHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
90578910aadSHong Zhang 
90678910aadSHong Zhang       while (i < k) {
90778910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
90878910aadSHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
90978910aadSHong Zhang         ili   = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
91078910aadSHong Zhang         uikdi = -ba[ili] * ba[bi[i]];
91178910aadSHong Zhang         dk += uikdi * ba[ili];
91278910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
91378910aadSHong Zhang 
91478910aadSHong Zhang         /* add multiple of row i to k-th row ... */
91578910aadSHong Zhang         jmin = ili + 1;
91678910aadSHong Zhang         nz   = bi[i + 1] - jmin;
91778910aadSHong Zhang         if (nz > 0) {
91878910aadSHong Zhang           bcol = bj + jmin;
91978910aadSHong Zhang           bval = ba + jmin;
92078910aadSHong Zhang           while (nz--) rtmp[*bcol++] += uikdi * (*bval++);
92178910aadSHong Zhang           /* update il and jl for i-th row */
92278910aadSHong Zhang           il[i] = jmin;
9239371c9d4SSatish Balay           j     = bj[jmin];
9249371c9d4SSatish Balay           jl[i] = jl[j];
9259371c9d4SSatish Balay           jl[j] = i;
92678910aadSHong Zhang         }
92778910aadSHong Zhang         i = nexti;
92878910aadSHong Zhang       }
92978910aadSHong Zhang 
9303cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
9313cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
9323cea4cbeSHong Zhang       rs   = 0.0;
93378910aadSHong Zhang       jmin = bi[k] + 1;
93478910aadSHong Zhang       nz   = bi[k + 1] - jmin;
93578910aadSHong Zhang       if (nz) {
93678910aadSHong Zhang         bcol = bj + jmin;
93778910aadSHong Zhang         while (nz--) {
93889f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
93989f845c8SHong Zhang           bcol++;
94078910aadSHong Zhang         }
94178910aadSHong Zhang       }
9423cea4cbeSHong Zhang 
9433cea4cbeSHong Zhang       sctx.rs = rs;
9443cea4cbeSHong Zhang       sctx.pv = dk;
9459566063dSJacob Faibussowitsch       PetscCall(MatPivotCheck(C, A, info, &sctx, k));
94607b50cabSHong Zhang       if (sctx.newshift) break; /* sctx.shift_amount is updated */
9470e95ead3SHong Zhang       dk = sctx.pv;
94878910aadSHong Zhang 
94978910aadSHong Zhang       /* copy data into U(k,:) */
95078910aadSHong Zhang       ba[bi[k]] = 1.0 / dk;
95178910aadSHong Zhang       jmin      = bi[k] + 1;
95278910aadSHong Zhang       nz        = bi[k + 1] - jmin;
95378910aadSHong Zhang       if (nz) {
95478910aadSHong Zhang         bcol = bj + jmin;
95578910aadSHong Zhang         bval = ba + jmin;
95678910aadSHong Zhang         while (nz--) {
95778910aadSHong Zhang           *bval++       = rtmp[*bcol];
95878910aadSHong Zhang           rtmp[*bcol++] = 0.0;
95978910aadSHong Zhang         }
96078910aadSHong Zhang         /* add k-th row into il and jl */
96178910aadSHong Zhang         il[k] = jmin;
9629371c9d4SSatish Balay         i     = bj[jmin];
9639371c9d4SSatish Balay         jl[k] = jl[i];
9649371c9d4SSatish Balay         jl[i] = k;
96578910aadSHong Zhang       }
96678910aadSHong Zhang     }
96707b50cabSHong Zhang   } while (sctx.newshift);
9689566063dSJacob Faibussowitsch   PetscCall(PetscFree3(rtmp, il, jl));
96978910aadSHong Zhang 
9700a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
9710a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
97278910aadSHong Zhang   C->assembled           = PETSC_TRUE;
97378910aadSHong Zhang   C->preallocated        = PETSC_TRUE;
97426fbe8dcSKarl Rupp 
9759566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(C->rmap->N));
9763cea4cbeSHong Zhang   if (sctx.nshift) {
977f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
9789566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shiftnz tries %" PetscInt_FMT ", shift_amount %g\n", sctx.nshift, (double)sctx.shift_amount));
979f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
9809566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shiftpd tries %" PetscInt_FMT ", shift_amount %g\n", sctx.nshift, (double)sctx.shift_amount));
98178910aadSHong Zhang     }
98278910aadSHong Zhang   }
983c05c3958SHong Zhang   PetscFunctionReturn(0);
984c05c3958SHong Zhang }
985c05c3958SHong Zhang 
986c6db04a5SJed Brown #include <petscbt.h>
987c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
9889371c9d4SSatish Balay PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact, Mat A, IS perm, const MatFactorInfo *info) {
98978910aadSHong Zhang   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data;
99078910aadSHong Zhang   Mat_SeqSBAIJ      *b;
99178910aadSHong Zhang   Mat                B;
99248dd3d27SHong Zhang   PetscBool          perm_identity, missing;
9935d0c19d7SBarry Smith   PetscInt           reallocs = 0, i, *ai = a->i, *aj = a->j, am = a->mbs, bs = A->rmap->bs, *ui;
9945d0c19d7SBarry Smith   const PetscInt    *rip;
99578910aadSHong Zhang   PetscInt           jmin, jmax, nzk, k, j, *jl, prow, *il, nextprow;
9960298fd71SBarry Smith   PetscInt           nlnk, *lnk, *lnk_lvl = NULL, ncols, ncols_upper, *cols, *cols_lvl, *uj, **uj_ptr, **uj_lvl_ptr;
99778910aadSHong Zhang   PetscReal          fill = info->fill, levels = info->levels;
9980298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL, current_space = NULL;
9990298fd71SBarry Smith   PetscFreeSpaceList free_space_lvl = NULL, current_space_lvl = NULL;
100078910aadSHong Zhang   PetscBT            lnkbt;
100178910aadSHong Zhang 
1002c05c3958SHong Zhang   PetscFunctionBegin;
10039566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal(A, &missing, &i));
100428b400f6SJacob Faibussowitsch   PetscCheck(!missing, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal entry %" PetscInt_FMT, i);
100548dd3d27SHong Zhang 
10066ad2eaddSHong Zhang   if (bs > 1) {
1007*48a46eb9SPierre Jolivet     if (!a->sbaijMat) PetscCall(MatConvert(A, MATSEQSBAIJ, MAT_INITIAL_MATRIX, &a->sbaijMat));
1008719d5645SBarry Smith     (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */
100926fbe8dcSKarl Rupp 
10109566063dSJacob Faibussowitsch     PetscCall(MatICCFactorSymbolic(fact, a->sbaijMat, perm, info));
10116ad2eaddSHong Zhang     PetscFunctionReturn(0);
10126ad2eaddSHong Zhang   }
10136ad2eaddSHong Zhang 
10149566063dSJacob Faibussowitsch   PetscCall(ISIdentity(perm, &perm_identity));
10159566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(perm, &rip));
101678910aadSHong Zhang 
101778910aadSHong Zhang   /* special case that simply copies fill pattern */
101878910aadSHong Zhang   if (!levels && perm_identity) {
10199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(am + 1, &ui));
102026fbe8dcSKarl Rupp     for (i = 0; i < am; i++) ui[i] = ai[i + 1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
1021719d5645SBarry Smith     B = fact;
10229566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(B, 1, 0, ui));
102378910aadSHong Zhang 
102478910aadSHong Zhang     b  = (Mat_SeqSBAIJ *)B->data;
102578910aadSHong Zhang     uj = b->j;
102678910aadSHong Zhang     for (i = 0; i < am; i++) {
102778910aadSHong Zhang       aj = a->j + a->diag[i];
102826fbe8dcSKarl Rupp       for (j = 0; j < ui[i]; j++) *uj++ = *aj++;
102978910aadSHong Zhang       b->ilen[i] = ui[i];
103078910aadSHong Zhang     }
10319566063dSJacob Faibussowitsch     PetscCall(PetscFree(ui));
103226fbe8dcSKarl Rupp 
1033d5f3da31SBarry Smith     B->factortype = MAT_FACTOR_NONE;
103426fbe8dcSKarl Rupp 
10359566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
10369566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
1037d5f3da31SBarry Smith     B->factortype = MAT_FACTOR_ICC;
103878910aadSHong Zhang 
103978910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
104078910aadSHong Zhang     PetscFunctionReturn(0);
104178910aadSHong Zhang   }
104278910aadSHong Zhang 
104378910aadSHong Zhang   /* initialization */
10449566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(am + 1, &ui));
1045e60cf9a0SBarry Smith   ui[0] = 0;
10469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * am + 1, &cols_lvl));
104778910aadSHong Zhang 
104878910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
104978910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
10509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(am, &uj_ptr, am, &uj_lvl_ptr, am, &il, am, &jl));
105178910aadSHong Zhang   for (i = 0; i < am; i++) {
10529371c9d4SSatish Balay     jl[i] = am;
10539371c9d4SSatish Balay     il[i] = 0;
105478910aadSHong Zhang   }
105578910aadSHong Zhang 
105678910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
105778910aadSHong Zhang   nlnk = am + 1;
10589566063dSJacob Faibussowitsch   PetscCall(PetscIncompleteLLCreate(am, am, nlnk, lnk, lnk_lvl, lnkbt));
105978910aadSHong Zhang 
106095b5ac22SHong Zhang   /* initial FreeSpace size is fill*(ai[am]+am)/2 */
10619566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill, PetscIntSumTruncate(ai[am] / 2, am / 2)), &free_space));
106226fbe8dcSKarl Rupp 
106378910aadSHong Zhang   current_space = free_space;
106426fbe8dcSKarl Rupp 
10659566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill, PetscIntSumTruncate(ai[am] / 2, am / 2)), &free_space_lvl));
106678910aadSHong Zhang   current_space_lvl = free_space_lvl;
106778910aadSHong Zhang 
106878910aadSHong Zhang   for (k = 0; k < am; k++) { /* for each active row k */
106978910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
107078910aadSHong Zhang     nzk         = 0;
107178910aadSHong Zhang     ncols       = ai[rip[k] + 1] - ai[rip[k]];
107278910aadSHong Zhang     ncols_upper = 0;
107378910aadSHong Zhang     cols        = cols_lvl + am;
107478910aadSHong Zhang     for (j = 0; j < ncols; j++) {
107578910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
107678910aadSHong Zhang       if (i >= k) { /* only take upper triangular entry */
107778910aadSHong Zhang         cols[ncols_upper]     = i;
107878910aadSHong Zhang         cols_lvl[ncols_upper] = -1; /* initialize level for nonzero entries */
107978910aadSHong Zhang         ncols_upper++;
108078910aadSHong Zhang       }
108178910aadSHong Zhang     }
10829566063dSJacob Faibussowitsch     PetscCall(PetscIncompleteLLAdd(ncols_upper, cols, levels, cols_lvl, am, &nlnk, lnk, lnk_lvl, lnkbt));
108378910aadSHong Zhang     nzk += nlnk;
108478910aadSHong Zhang 
108578910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
108678910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
108778910aadSHong Zhang 
108878910aadSHong Zhang     while (prow < k) {
108978910aadSHong Zhang       nextprow = jl[prow];
109078910aadSHong Zhang 
109178910aadSHong Zhang       /* merge prow into k-th row */
109278910aadSHong Zhang       jmin  = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */
109378910aadSHong Zhang       jmax  = ui[prow + 1];
109478910aadSHong Zhang       ncols = jmax - jmin;
109578910aadSHong Zhang       i     = jmin - ui[prow];
109678910aadSHong Zhang       cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
109778910aadSHong Zhang       for (j = 0; j < ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
10989566063dSJacob Faibussowitsch       PetscCall(PetscIncompleteLLAddSorted(ncols, cols, levels, cols_lvl, am, &nlnk, lnk, lnk_lvl, lnkbt));
109978910aadSHong Zhang       nzk += nlnk;
110078910aadSHong Zhang 
110178910aadSHong Zhang       /* update il and jl for prow */
110278910aadSHong Zhang       if (jmin < jmax) {
110378910aadSHong Zhang         il[prow] = jmin;
110426fbe8dcSKarl Rupp 
11059371c9d4SSatish Balay         j        = *cols;
11069371c9d4SSatish Balay         jl[prow] = jl[j];
11079371c9d4SSatish Balay         jl[j]    = prow;
110878910aadSHong Zhang       }
110978910aadSHong Zhang       prow = nextprow;
111078910aadSHong Zhang     }
111178910aadSHong Zhang 
111278910aadSHong Zhang     /* if free space is not available, make more free space */
111378910aadSHong Zhang     if (current_space->local_remaining < nzk) {
111478910aadSHong Zhang       i = am - k + 1;                                                             /* num of unfactored rows */
1115f91af8c7SBarry Smith       i = PetscMin(PetscIntMultTruncate(i, nzk), PetscIntMultTruncate(i, i - 1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
11169566063dSJacob Faibussowitsch       PetscCall(PetscFreeSpaceGet(i, &current_space));
11179566063dSJacob Faibussowitsch       PetscCall(PetscFreeSpaceGet(i, &current_space_lvl));
111878910aadSHong Zhang       reallocs++;
111978910aadSHong Zhang     }
112078910aadSHong Zhang 
112178910aadSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
11229566063dSJacob Faibussowitsch     PetscCall(PetscIncompleteLLClean(am, am, nzk, lnk, lnk_lvl, current_space->array, current_space_lvl->array, lnkbt));
112378910aadSHong Zhang 
112478910aadSHong Zhang     /* add the k-th row into il and jl */
112578910aadSHong Zhang     if (nzk - 1 > 0) {
112678910aadSHong Zhang       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
11279371c9d4SSatish Balay       jl[k] = jl[i];
11289371c9d4SSatish Balay       jl[i] = k;
112978910aadSHong Zhang       il[k] = ui[k] + 1;
113078910aadSHong Zhang     }
113178910aadSHong Zhang     uj_ptr[k]     = current_space->array;
113278910aadSHong Zhang     uj_lvl_ptr[k] = current_space_lvl->array;
113378910aadSHong Zhang 
113478910aadSHong Zhang     current_space->array += nzk;
113578910aadSHong Zhang     current_space->local_used += nzk;
113678910aadSHong Zhang     current_space->local_remaining -= nzk;
113778910aadSHong Zhang 
113878910aadSHong Zhang     current_space_lvl->array += nzk;
113978910aadSHong Zhang     current_space_lvl->local_used += nzk;
114078910aadSHong Zhang     current_space_lvl->local_remaining -= nzk;
114178910aadSHong Zhang 
114278910aadSHong Zhang     ui[k + 1] = ui[k] + nzk;
114378910aadSHong Zhang   }
114478910aadSHong Zhang 
11459566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(perm, &rip));
11469566063dSJacob Faibussowitsch   PetscCall(PetscFree4(uj_ptr, uj_lvl_ptr, il, jl));
11479566063dSJacob Faibussowitsch   PetscCall(PetscFree(cols_lvl));
114878910aadSHong Zhang 
11499263d837SHong Zhang   /* copy free_space into uj and free free_space; set uj in new datastructure; */
11509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ui[am] + 1, &uj));
11519566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceContiguous(&free_space, uj));
11529566063dSJacob Faibussowitsch   PetscCall(PetscIncompleteLLDestroy(lnk, lnkbt));
11539566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceDestroy(free_space_lvl));
115478910aadSHong Zhang 
115578910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1156719d5645SBarry Smith   B = fact;
11579566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(B, 1, MAT_SKIP_ALLOCATION, NULL));
115878910aadSHong Zhang 
115978910aadSHong Zhang   b               = (Mat_SeqSBAIJ *)B->data;
116078910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
1161e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1162e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
116326fbe8dcSKarl Rupp 
11649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ui[am] + 1, &b->a));
116526fbe8dcSKarl Rupp 
116678910aadSHong Zhang   b->j             = uj;
116778910aadSHong Zhang   b->i             = ui;
1168f4259b30SLisandro Dalcin   b->diag          = NULL;
1169f4259b30SLisandro Dalcin   b->ilen          = NULL;
1170f4259b30SLisandro Dalcin   b->imax          = NULL;
117178910aadSHong Zhang   b->row           = perm;
117278910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
117326fbe8dcSKarl Rupp 
11749566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)perm));
117526fbe8dcSKarl Rupp 
117678910aadSHong Zhang   b->icol = perm;
117726fbe8dcSKarl Rupp 
11789566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)perm));
11799566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(am + 1, &b->solve_work));
11809566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)B, (ui[am] - am) * (sizeof(PetscInt) + sizeof(MatScalar))));
118126fbe8dcSKarl Rupp 
118278910aadSHong Zhang   b->maxnz = b->nz = ui[am];
118378910aadSHong Zhang 
118478910aadSHong Zhang   B->info.factor_mallocs   = reallocs;
118578910aadSHong Zhang   B->info.fill_ratio_given = fill;
118675567043SBarry Smith   if (ai[am] != 0.) {
118795b5ac22SHong Zhang     /* nonzeros in lower triangular part of A (includign diagonals)= (ai[am]+am)/2 */
118895b5ac22SHong Zhang     B->info.fill_ratio_needed = ((PetscReal)2 * ui[am]) / (ai[am] + am);
118978910aadSHong Zhang   } else {
119078910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
119178910aadSHong Zhang   }
11929263d837SHong Zhang #if defined(PETSC_USE_INFO)
11939263d837SHong Zhang   if (ai[am] != 0) {
119495b5ac22SHong Zhang     PetscReal af = B->info.fill_ratio_needed;
11959566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n", reallocs, (double)fill, (double)af));
11969566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Run with -pc_factor_fill %g or use \n", (double)af));
11979566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "PCFactorSetFill(pc,%g) for best performance.\n", (double)af));
11989263d837SHong Zhang   } else {
11999566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Empty matrix\n"));
12009263d837SHong Zhang   }
12019263d837SHong Zhang #endif
120278910aadSHong Zhang   if (perm_identity) {
12030a3c351aSHong Zhang     B->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
12040a3c351aSHong Zhang     B->ops->solvetranspose        = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
120578910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
120678910aadSHong Zhang   } else {
1207719d5645SBarry Smith     (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
120878910aadSHong Zhang   }
1209c05c3958SHong Zhang   PetscFunctionReturn(0);
1210c05c3958SHong Zhang }
1211c05c3958SHong Zhang 
12129371c9d4SSatish Balay PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact, Mat A, IS perm, const MatFactorInfo *info) {
121378910aadSHong Zhang   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data;
121478910aadSHong Zhang   Mat_SeqSBAIJ      *b;
121578910aadSHong Zhang   Mat                B;
12169186b105SHong Zhang   PetscBool          perm_identity, missing;
121778910aadSHong Zhang   PetscReal          fill = info->fill;
12185d0c19d7SBarry Smith   const PetscInt    *rip;
12195d0c19d7SBarry Smith   PetscInt           i, mbs = a->mbs, bs = A->rmap->bs, *ai = a->i, *aj = a->j, reallocs = 0, prow;
122078910aadSHong Zhang   PetscInt          *jl, jmin, jmax, nzk, *ui, k, j, *il, nextprow;
122178910aadSHong Zhang   PetscInt           nlnk, *lnk, ncols, ncols_upper, *cols, *uj, **ui_ptr, *uj_ptr;
12220298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL, current_space = NULL;
122378910aadSHong Zhang   PetscBT            lnkbt;
122478910aadSHong Zhang 
1225c05c3958SHong Zhang   PetscFunctionBegin;
12266ad2eaddSHong Zhang   if (bs > 1) { /* convert to seqsbaij */
1227*48a46eb9SPierre Jolivet     if (!a->sbaijMat) PetscCall(MatConvert(A, MATSEQSBAIJ, MAT_INITIAL_MATRIX, &a->sbaijMat));
1228719d5645SBarry Smith     (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */
122926fbe8dcSKarl Rupp 
12309566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorSymbolic(fact, a->sbaijMat, perm, info));
12316ad2eaddSHong Zhang     PetscFunctionReturn(0);
12326ad2eaddSHong Zhang   }
12336ad2eaddSHong Zhang 
12349566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal(A, &missing, &i));
123528b400f6SJacob Faibussowitsch   PetscCheck(!missing, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal entry %" PetscInt_FMT, i);
12369186b105SHong Zhang 
123778910aadSHong Zhang   /* check whether perm is the identity mapping */
12389566063dSJacob Faibussowitsch   PetscCall(ISIdentity(perm, &perm_identity));
123928b400f6SJacob Faibussowitsch   PetscCheck(perm_identity, PETSC_COMM_SELF, PETSC_ERR_SUP, "Matrix reordering is not supported");
12409566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(perm, &rip));
124178910aadSHong Zhang 
124278910aadSHong Zhang   /* initialization */
12439566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs + 1, &ui));
1244e60cf9a0SBarry Smith   ui[0] = 0;
124578910aadSHong Zhang 
124678910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
124778910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
12489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(mbs, &ui_ptr, mbs, &il, mbs, &jl, mbs, &cols));
124978910aadSHong Zhang   for (i = 0; i < mbs; i++) {
12509371c9d4SSatish Balay     jl[i] = mbs;
12519371c9d4SSatish Balay     il[i] = 0;
125278910aadSHong Zhang   }
125378910aadSHong Zhang 
125478910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
125578910aadSHong Zhang   nlnk = mbs + 1;
12569566063dSJacob Faibussowitsch   PetscCall(PetscLLCreate(mbs, mbs, nlnk, lnk, lnkbt));
125778910aadSHong Zhang 
12586a69fef8SHong Zhang   /* initial FreeSpace size is fill* (ai[mbs]+mbs)/2 */
12599566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill, PetscIntSumTruncate(ai[mbs] / 2, mbs / 2)), &free_space));
126026fbe8dcSKarl Rupp 
126178910aadSHong Zhang   current_space = free_space;
126278910aadSHong Zhang 
126378910aadSHong Zhang   for (k = 0; k < mbs; k++) { /* for each active row k */
126478910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
126578910aadSHong Zhang     nzk         = 0;
126678910aadSHong Zhang     ncols       = ai[rip[k] + 1] - ai[rip[k]];
1267e60cf9a0SBarry Smith     ncols_upper = 0;
126878910aadSHong Zhang     for (j = 0; j < ncols; j++) {
126978910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
127078910aadSHong Zhang       if (i >= k) { /* only take upper triangular entry */
127178910aadSHong Zhang         cols[ncols_upper] = i;
127278910aadSHong Zhang         ncols_upper++;
127378910aadSHong Zhang       }
127478910aadSHong Zhang     }
12759566063dSJacob Faibussowitsch     PetscCall(PetscLLAdd(ncols_upper, cols, mbs, &nlnk, lnk, lnkbt));
127678910aadSHong Zhang     nzk += nlnk;
127778910aadSHong Zhang 
127878910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
127978910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
128078910aadSHong Zhang 
128178910aadSHong Zhang     while (prow < k) {
128278910aadSHong Zhang       nextprow = jl[prow];
128378910aadSHong Zhang       /* merge prow into k-th row */
128478910aadSHong Zhang       jmin     = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
128578910aadSHong Zhang       jmax     = ui[prow + 1];
128678910aadSHong Zhang       ncols    = jmax - jmin;
128778910aadSHong Zhang       uj_ptr   = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
12889566063dSJacob Faibussowitsch       PetscCall(PetscLLAddSorted(ncols, uj_ptr, mbs, &nlnk, lnk, lnkbt));
128978910aadSHong Zhang       nzk += nlnk;
129078910aadSHong Zhang 
129178910aadSHong Zhang       /* update il and jl for prow */
129278910aadSHong Zhang       if (jmin < jmax) {
129378910aadSHong Zhang         il[prow] = jmin;
129426fbe8dcSKarl Rupp         j        = *uj_ptr;
129526fbe8dcSKarl Rupp         jl[prow] = jl[j];
129626fbe8dcSKarl Rupp         jl[j]    = prow;
129778910aadSHong Zhang       }
129878910aadSHong Zhang       prow = nextprow;
129978910aadSHong Zhang     }
130078910aadSHong Zhang 
130178910aadSHong Zhang     /* if free space is not available, make more free space */
130278910aadSHong Zhang     if (current_space->local_remaining < nzk) {
130378910aadSHong Zhang       i = mbs - k + 1;                                                            /* num of unfactored rows */
1304f91af8c7SBarry Smith       i = PetscMin(PetscIntMultTruncate(i, nzk), PetscIntMultTruncate(i, i - 1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
13059566063dSJacob Faibussowitsch       PetscCall(PetscFreeSpaceGet(i, &current_space));
130678910aadSHong Zhang       reallocs++;
130778910aadSHong Zhang     }
130878910aadSHong Zhang 
130978910aadSHong Zhang     /* copy data into free space, then initialize lnk */
13109566063dSJacob Faibussowitsch     PetscCall(PetscLLClean(mbs, mbs, nzk, lnk, current_space->array, lnkbt));
131178910aadSHong Zhang 
131278910aadSHong Zhang     /* add the k-th row into il and jl */
131378910aadSHong Zhang     if (nzk - 1 > 0) {
131478910aadSHong Zhang       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
13159371c9d4SSatish Balay       jl[k] = jl[i];
13169371c9d4SSatish Balay       jl[i] = k;
131778910aadSHong Zhang       il[k] = ui[k] + 1;
131878910aadSHong Zhang     }
131978910aadSHong Zhang     ui_ptr[k] = current_space->array;
132078910aadSHong Zhang     current_space->array += nzk;
132178910aadSHong Zhang     current_space->local_used += nzk;
132278910aadSHong Zhang     current_space->local_remaining -= nzk;
132378910aadSHong Zhang 
132478910aadSHong Zhang     ui[k + 1] = ui[k] + nzk;
132578910aadSHong Zhang   }
132678910aadSHong Zhang 
13279566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(perm, &rip));
13289566063dSJacob Faibussowitsch   PetscCall(PetscFree4(ui_ptr, il, jl, cols));
132978910aadSHong Zhang 
13309263d837SHong Zhang   /* copy free_space into uj and free free_space; set uj in new datastructure; */
13319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ui[mbs] + 1, &uj));
13329566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceContiguous(&free_space, uj));
13339566063dSJacob Faibussowitsch   PetscCall(PetscLLDestroy(lnk, lnkbt));
133478910aadSHong Zhang 
133578910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1336719d5645SBarry Smith   B = fact;
13379566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(B, bs, MAT_SKIP_ALLOCATION, NULL));
133878910aadSHong Zhang 
133978910aadSHong Zhang   b               = (Mat_SeqSBAIJ *)B->data;
134078910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
1341e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1342e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
134326fbe8dcSKarl Rupp 
13449566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ui[mbs] + 1, &b->a));
134526fbe8dcSKarl Rupp 
134678910aadSHong Zhang   b->j             = uj;
134778910aadSHong Zhang   b->i             = ui;
1348f4259b30SLisandro Dalcin   b->diag          = NULL;
1349f4259b30SLisandro Dalcin   b->ilen          = NULL;
1350f4259b30SLisandro Dalcin   b->imax          = NULL;
135178910aadSHong Zhang   b->row           = perm;
135278910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
135326fbe8dcSKarl Rupp 
13549566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)perm));
135578910aadSHong Zhang   b->icol = perm;
13569566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)perm));
13579566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs + 1, &b->solve_work));
13589566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)B, (ui[mbs] - mbs) * (sizeof(PetscInt) + sizeof(MatScalar))));
135978910aadSHong Zhang   b->maxnz = b->nz = ui[mbs];
136078910aadSHong Zhang 
136178910aadSHong Zhang   B->info.factor_mallocs   = reallocs;
136278910aadSHong Zhang   B->info.fill_ratio_given = fill;
136375567043SBarry Smith   if (ai[mbs] != 0.) {
136495b5ac22SHong Zhang     /* nonzeros in lower triangular part of A = (ai[mbs]+mbs)/2 */
136595b5ac22SHong Zhang     B->info.fill_ratio_needed = ((PetscReal)2 * ui[mbs]) / (ai[mbs] + mbs);
136678910aadSHong Zhang   } else {
136778910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
136878910aadSHong Zhang   }
13699263d837SHong Zhang #if defined(PETSC_USE_INFO)
13709263d837SHong Zhang   if (ai[mbs] != 0.) {
13719263d837SHong Zhang     PetscReal af = B->info.fill_ratio_needed;
13729566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n", reallocs, (double)fill, (double)af));
13739566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Run with -pc_factor_fill %g or use \n", (double)af));
13749566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "PCFactorSetFill(pc,%g) for best performance.\n", (double)af));
13759263d837SHong Zhang   } else {
13769566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Empty matrix\n"));
13779263d837SHong Zhang   }
13789263d837SHong Zhang #endif
137978910aadSHong Zhang   if (perm_identity) {
13806ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
138178910aadSHong Zhang   } else {
13826ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
138378910aadSHong Zhang   }
1384c05c3958SHong Zhang   PetscFunctionReturn(0);
1385c05c3958SHong Zhang }
1386c8342467SHong Zhang 
13879371c9d4SSatish Balay PetscErrorCode MatSolve_SeqBAIJ_N_NaturalOrdering(Mat A, Vec bb, Vec xx) {
13881a83e813SShri Abhyankar   Mat_SeqBAIJ       *a  = (Mat_SeqBAIJ *)A->data;
13891a83e813SShri Abhyankar   const PetscInt    *ai = a->i, *aj = a->j, *adiag = a->diag, *vi;
13901a83e813SShri Abhyankar   PetscInt           i, k, n                       = a->mbs;
13911a83e813SShri Abhyankar   PetscInt           nz, bs = A->rmap->bs, bs2 = a->bs2;
1392d9ca1df4SBarry Smith   const MatScalar   *aa = a->a, *v;
1393d9ca1df4SBarry Smith   PetscScalar       *x, *s, *t, *ls;
1394d9ca1df4SBarry Smith   const PetscScalar *b;
13951a83e813SShri Abhyankar 
13961a83e813SShri Abhyankar   PetscFunctionBegin;
13979566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
13989566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
13991a83e813SShri Abhyankar   t = a->solve_work;
14001a83e813SShri Abhyankar 
14011a83e813SShri Abhyankar   /* forward solve the lower triangular */
14029566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t, b, bs)); /* copy 1st block of b to t */
14031a83e813SShri Abhyankar 
14041a83e813SShri Abhyankar   for (i = 1; i < n; i++) {
14051a83e813SShri Abhyankar     v  = aa + bs2 * ai[i];
14061a83e813SShri Abhyankar     vi = aj + ai[i];
14071a83e813SShri Abhyankar     nz = ai[i + 1] - ai[i];
14081a83e813SShri Abhyankar     s  = t + bs * i;
14099566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(s, b + bs * i, bs)); /* copy i_th block of b to t */
14101a83e813SShri Abhyankar     for (k = 0; k < nz; k++) {
141196b95a6bSBarry Smith       PetscKernel_v_gets_v_minus_A_times_w(bs, s, v, t + bs * vi[k]);
14121a83e813SShri Abhyankar       v += bs2;
14131a83e813SShri Abhyankar     }
14141a83e813SShri Abhyankar   }
14151a83e813SShri Abhyankar 
14161a83e813SShri Abhyankar   /* backward solve the upper triangular */
14171a83e813SShri Abhyankar   ls = a->solve_work + A->cmap->n;
14181a83e813SShri Abhyankar   for (i = n - 1; i >= 0; i--) {
14191a83e813SShri Abhyankar     v  = aa + bs2 * (adiag[i + 1] + 1);
14201a83e813SShri Abhyankar     vi = aj + adiag[i + 1] + 1;
14211a83e813SShri Abhyankar     nz = adiag[i] - adiag[i + 1] - 1;
14229566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(ls, t + i * bs, bs));
14231a83e813SShri Abhyankar     for (k = 0; k < nz; k++) {
142496b95a6bSBarry Smith       PetscKernel_v_gets_v_minus_A_times_w(bs, ls, v, t + bs * vi[k]);
14251a83e813SShri Abhyankar       v += bs2;
14261a83e813SShri Abhyankar     }
142796b95a6bSBarry Smith     PetscKernel_w_gets_A_times_v(bs, ls, aa + bs2 * adiag[i], t + i * bs); /* *inv(diagonal[i]) */
14289566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(x + i * bs, t + i * bs, bs));
14291a83e813SShri Abhyankar   }
14301a83e813SShri Abhyankar 
14319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
14329566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
14339566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * (a->bs2) * (a->nz) - A->rmap->bs * A->cmap->n));
14341a83e813SShri Abhyankar   PetscFunctionReturn(0);
14351a83e813SShri Abhyankar }
14361a83e813SShri Abhyankar 
14379371c9d4SSatish Balay PetscErrorCode MatSolve_SeqBAIJ_N(Mat A, Vec bb, Vec xx) {
143835aa4fcfSShri Abhyankar   Mat_SeqBAIJ       *a     = (Mat_SeqBAIJ *)A->data;
143935aa4fcfSShri Abhyankar   IS                 iscol = a->col, isrow = a->row;
144035aa4fcfSShri Abhyankar   const PetscInt    *r, *c, *rout, *cout, *ai = a->i, *aj = a->j, *adiag = a->diag, *vi;
144135aa4fcfSShri Abhyankar   PetscInt           i, m, n = a->mbs;
144235aa4fcfSShri Abhyankar   PetscInt           nz, bs = A->rmap->bs, bs2 = a->bs2;
1443d9ca1df4SBarry Smith   const MatScalar   *aa = a->a, *v;
1444d9ca1df4SBarry Smith   PetscScalar       *x, *s, *t, *ls;
1445d9ca1df4SBarry Smith   const PetscScalar *b;
144635aa4fcfSShri Abhyankar 
144735aa4fcfSShri Abhyankar   PetscFunctionBegin;
14489566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
14499566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
145035aa4fcfSShri Abhyankar   t = a->solve_work;
145135aa4fcfSShri Abhyankar 
14529371c9d4SSatish Balay   PetscCall(ISGetIndices(isrow, &rout));
14539371c9d4SSatish Balay   r = rout;
14549371c9d4SSatish Balay   PetscCall(ISGetIndices(iscol, &cout));
14559371c9d4SSatish Balay   c = cout;
145635aa4fcfSShri Abhyankar 
145735aa4fcfSShri Abhyankar   /* forward solve the lower triangular */
14589566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t, b + bs * r[0], bs));
145935aa4fcfSShri Abhyankar   for (i = 1; i < n; i++) {
146035aa4fcfSShri Abhyankar     v  = aa + bs2 * ai[i];
146135aa4fcfSShri Abhyankar     vi = aj + ai[i];
146235aa4fcfSShri Abhyankar     nz = ai[i + 1] - ai[i];
146335aa4fcfSShri Abhyankar     s  = t + bs * i;
14649566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(s, b + bs * r[i], bs));
146535aa4fcfSShri Abhyankar     for (m = 0; m < nz; m++) {
146696b95a6bSBarry Smith       PetscKernel_v_gets_v_minus_A_times_w(bs, s, v, t + bs * vi[m]);
146735aa4fcfSShri Abhyankar       v += bs2;
146835aa4fcfSShri Abhyankar     }
146935aa4fcfSShri Abhyankar   }
147035aa4fcfSShri Abhyankar 
147135aa4fcfSShri Abhyankar   /* backward solve the upper triangular */
147235aa4fcfSShri Abhyankar   ls = a->solve_work + A->cmap->n;
147335aa4fcfSShri Abhyankar   for (i = n - 1; i >= 0; i--) {
147435aa4fcfSShri Abhyankar     v  = aa + bs2 * (adiag[i + 1] + 1);
147535aa4fcfSShri Abhyankar     vi = aj + adiag[i + 1] + 1;
147635aa4fcfSShri Abhyankar     nz = adiag[i] - adiag[i + 1] - 1;
14779566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(ls, t + i * bs, bs));
147835aa4fcfSShri Abhyankar     for (m = 0; m < nz; m++) {
147996b95a6bSBarry Smith       PetscKernel_v_gets_v_minus_A_times_w(bs, ls, v, t + bs * vi[m]);
148035aa4fcfSShri Abhyankar       v += bs2;
148135aa4fcfSShri Abhyankar     }
148296b95a6bSBarry Smith     PetscKernel_w_gets_A_times_v(bs, ls, v, t + i * bs); /* *inv(diagonal[i]) */
14839566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(x + bs * c[i], t + i * bs, bs));
148435aa4fcfSShri Abhyankar   }
14859566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &rout));
14869566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscol, &cout));
14879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
14889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
14899566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * (a->bs2) * (a->nz) - A->rmap->bs * A->cmap->n));
149035aa4fcfSShri Abhyankar   PetscFunctionReturn(0);
149135aa4fcfSShri Abhyankar }
149235aa4fcfSShri Abhyankar 
1493a25532f0SBarry Smith /*
1494a25532f0SBarry Smith     For each block in an block array saves the largest absolute value in the block into another array
1495a25532f0SBarry Smith */
14969371c9d4SSatish Balay static PetscErrorCode MatBlockAbs_private(PetscInt nbs, PetscInt bs2, PetscScalar *blockarray, PetscReal *absarray) {
14972efa7f71SHong Zhang   PetscInt i, j;
14985fd66863SKarl Rupp 
14992efa7f71SHong Zhang   PetscFunctionBegin;
15009566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(absarray, nbs + 1));
15012efa7f71SHong Zhang   for (i = 0; i < nbs; i++) {
15022efa7f71SHong Zhang     for (j = 0; j < bs2; j++) {
15032efa7f71SHong Zhang       if (absarray[i] < PetscAbsScalar(blockarray[i * nbs + j])) absarray[i] = PetscAbsScalar(blockarray[i * nbs + j]);
15042efa7f71SHong Zhang     }
15052efa7f71SHong Zhang   }
15062efa7f71SHong Zhang   PetscFunctionReturn(0);
15072efa7f71SHong Zhang }
15082efa7f71SHong Zhang 
1509fe97e370SBarry Smith /*
1510fe97e370SBarry Smith      This needs to be renamed and called by the regular MatILUFactor_SeqBAIJ when drop tolerance is used
1511fe97e370SBarry Smith */
15129371c9d4SSatish Balay PetscErrorCode MatILUDTFactor_SeqBAIJ(Mat A, IS isrow, IS iscol, const MatFactorInfo *info, Mat *fact) {
15132efa7f71SHong Zhang   Mat             B = *fact;
15142efa7f71SHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b;
15152efa7f71SHong Zhang   IS              isicol;
15162efa7f71SHong Zhang   const PetscInt *r, *ic;
15172efa7f71SHong Zhang   PetscInt        i, mbs = a->mbs, bs = A->rmap->bs, bs2 = a->bs2, *ai = a->i, *aj = a->j, *ajtmp, *adiag;
15182efa7f71SHong Zhang   PetscInt       *bi, *bj, *bdiag;
15192efa7f71SHong Zhang 
15202efa7f71SHong Zhang   PetscInt   row, nzi, nzi_bl, nzi_bu, *im, dtcount, nzi_al, nzi_au;
15212efa7f71SHong Zhang   PetscInt   nlnk, *lnk;
15222efa7f71SHong Zhang   PetscBT    lnkbt;
1523ace3abfcSBarry Smith   PetscBool  row_identity, icol_identity;
15242efa7f71SHong Zhang   MatScalar *aatmp, *pv, *batmp, *ba, *rtmp, *pc, *multiplier, *vtmp;
15252efa7f71SHong Zhang   PetscInt   j, nz, *pj, *bjtmp, k, ncut, *jtmp;
15262efa7f71SHong Zhang 
1527182b8fbaSHong Zhang   PetscReal  dt = info->dt; /* shift=info->shiftamount; */
15282efa7f71SHong Zhang   PetscInt   nnz_max;
1529ace3abfcSBarry Smith   PetscBool  missing;
1530021822bcSHong Zhang   PetscReal *vtmp_abs;
153197ef94ebSSatish Balay   MatScalar *v_work;
153297ef94ebSSatish Balay   PetscInt  *v_pivots;
15335f8bbccaSHong Zhang   PetscBool  allowzeropivot, zeropivotdetected = PETSC_FALSE;
15342efa7f71SHong Zhang 
1535c8342467SHong Zhang   PetscFunctionBegin;
15362efa7f71SHong Zhang   /* ------- symbolic factorization, can be reused ---------*/
15379566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal(A, &missing, &i));
153828b400f6SJacob Faibussowitsch   PetscCheck(!missing, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal entry %" PetscInt_FMT, i);
15392efa7f71SHong Zhang   adiag = a->diag;
15402efa7f71SHong Zhang 
15419566063dSJacob Faibussowitsch   PetscCall(ISInvertPermutation(iscol, PETSC_DECIDE, &isicol));
15422efa7f71SHong Zhang 
15432efa7f71SHong Zhang   /* bdiag is location of diagonal in factor */
15449566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs + 1, &bdiag));
15452efa7f71SHong Zhang 
15462efa7f71SHong Zhang   /* allocate row pointers bi */
15479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * mbs + 2, &bi));
15482efa7f71SHong Zhang 
15492efa7f71SHong Zhang   /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */
15502efa7f71SHong Zhang   dtcount = (PetscInt)info->dtcount;
15516bce7ff8SHong Zhang   if (dtcount > mbs - 1) dtcount = mbs - 1;
15526bce7ff8SHong Zhang   nnz_max = ai[mbs] + 2 * mbs * dtcount + 2;
15536da515a1SHong Zhang   /* printf("MatILUDTFactor_SeqBAIJ, bs %d, ai[mbs] %d, nnz_max  %d, dtcount %d\n",bs,ai[mbs],nnz_max,dtcount); */
15549566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nnz_max, &bj));
15556bce7ff8SHong Zhang   nnz_max = nnz_max * bs2;
15569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nnz_max, &ba));
15572efa7f71SHong Zhang 
15582efa7f71SHong Zhang   /* put together the new matrix */
15599566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(B, bs, MAT_SKIP_ALLOCATION, NULL));
15609566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)B, (PetscObject)isicol));
156126fbe8dcSKarl Rupp 
15622efa7f71SHong Zhang   b               = (Mat_SeqBAIJ *)(B)->data;
15632efa7f71SHong Zhang   b->free_a       = PETSC_TRUE;
15642efa7f71SHong Zhang   b->free_ij      = PETSC_TRUE;
15652efa7f71SHong Zhang   b->singlemalloc = PETSC_FALSE;
156626fbe8dcSKarl Rupp 
15672efa7f71SHong Zhang   b->a    = ba;
15682efa7f71SHong Zhang   b->j    = bj;
15692efa7f71SHong Zhang   b->i    = bi;
15702efa7f71SHong Zhang   b->diag = bdiag;
1571f4259b30SLisandro Dalcin   b->ilen = NULL;
1572f4259b30SLisandro Dalcin   b->imax = NULL;
15732efa7f71SHong Zhang   b->row  = isrow;
15742efa7f71SHong Zhang   b->col  = iscol;
157526fbe8dcSKarl Rupp 
15769566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)isrow));
15779566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)iscol));
157826fbe8dcSKarl Rupp 
15792efa7f71SHong Zhang   b->icol = isicol;
15809566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs * (mbs + 1), &b->solve_work));
15819566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)B, nnz_max * (sizeof(PetscInt) + sizeof(MatScalar))));
15822efa7f71SHong Zhang   b->maxnz = nnz_max / bs2;
15832efa7f71SHong Zhang 
1584d5f3da31SBarry Smith   (B)->factortype            = MAT_FACTOR_ILUDT;
15852efa7f71SHong Zhang   (B)->info.factor_mallocs   = 0;
15862efa7f71SHong Zhang   (B)->info.fill_ratio_given = ((PetscReal)nnz_max) / ((PetscReal)(ai[mbs] * bs2));
15872efa7f71SHong Zhang   /* ------- end of symbolic factorization ---------*/
15889566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
15899566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
15902efa7f71SHong Zhang 
15912efa7f71SHong Zhang   /* linked list for storing column indices of the active row */
15922efa7f71SHong Zhang   nlnk = mbs + 1;
15939566063dSJacob Faibussowitsch   PetscCall(PetscLLCreate(mbs, mbs, nlnk, lnk, lnkbt));
15942efa7f71SHong Zhang 
15952efa7f71SHong Zhang   /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */
15969566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(mbs, &im, mbs, &jtmp));
15972efa7f71SHong Zhang   /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */
15989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(mbs * bs2, &rtmp, mbs * bs2, &vtmp));
15999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs + 1, &vtmp_abs));
16009566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(bs, &v_work, bs2, &multiplier, bs, &v_pivots));
16012efa7f71SHong Zhang 
16025f8bbccaSHong Zhang   allowzeropivot  = PetscNot(A->erroriffailure);
1603e60cf9a0SBarry Smith   bi[0]           = 0;
16042efa7f71SHong Zhang   bdiag[0]        = (nnz_max / bs2) - 1; /* location of diagonal in factor B */
16056bce7ff8SHong Zhang   bi[2 * mbs + 1] = bdiag[0] + 1;        /* endof bj and ba array */
16062efa7f71SHong Zhang   for (i = 0; i < mbs; i++) {
16072efa7f71SHong Zhang     /* copy initial fill into linked list */
16082efa7f71SHong Zhang     nzi = ai[r[i] + 1] - ai[r[i]];
160928b400f6SJacob Faibussowitsch     PetscCheck(nzi, PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Empty row in matrix: row in original ordering %" PetscInt_FMT " in permuted ordering %" PetscInt_FMT, r[i], i);
16102efa7f71SHong Zhang     nzi_al = adiag[r[i]] - ai[r[i]];
16112efa7f71SHong Zhang     nzi_au = ai[r[i] + 1] - adiag[r[i]] - 1;
16122efa7f71SHong Zhang 
16132efa7f71SHong Zhang     /* load in initial unfactored row */
16142efa7f71SHong Zhang     ajtmp = aj + ai[r[i]];
16159566063dSJacob Faibussowitsch     PetscCall(PetscLLAddPerm(nzi, ajtmp, ic, mbs, &nlnk, lnk, lnkbt));
16169566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(rtmp, mbs * bs2));
16172efa7f71SHong Zhang     aatmp = a->a + bs2 * ai[r[i]];
16189566063dSJacob Faibussowitsch     for (j = 0; j < nzi; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ic[ajtmp[j]], aatmp + bs2 * j, bs2));
16192efa7f71SHong Zhang 
16202efa7f71SHong Zhang     /* add pivot rows into linked list */
16212efa7f71SHong Zhang     row = lnk[mbs];
16222efa7f71SHong Zhang     while (row < i) {
16232efa7f71SHong Zhang       nzi_bl = bi[row + 1] - bi[row] + 1;
16242efa7f71SHong Zhang       bjtmp  = bj + bdiag[row + 1] + 1; /* points to 1st column next to the diagonal in U */
16259566063dSJacob Faibussowitsch       PetscCall(PetscLLAddSortedLU(bjtmp, row, &nlnk, lnk, lnkbt, i, nzi_bl, im));
16262efa7f71SHong Zhang       nzi += nlnk;
16272efa7f71SHong Zhang       row = lnk[row];
16282efa7f71SHong Zhang     }
16292efa7f71SHong Zhang 
16302efa7f71SHong Zhang     /* copy data from lnk into jtmp, then initialize lnk */
16319566063dSJacob Faibussowitsch     PetscCall(PetscLLClean(mbs, mbs, nzi, lnk, jtmp, lnkbt));
16322efa7f71SHong Zhang 
16332efa7f71SHong Zhang     /* numerical factorization */
16342efa7f71SHong Zhang     bjtmp = jtmp;
16352efa7f71SHong Zhang     row   = *bjtmp++; /* 1st pivot row */
16362efa7f71SHong Zhang 
16372efa7f71SHong Zhang     while (row < i) {
16382efa7f71SHong Zhang       pc = rtmp + bs2 * row;
16392efa7f71SHong Zhang       pv = ba + bs2 * bdiag[row];                           /* inv(diag) of the pivot row */
164096b95a6bSBarry Smith       PetscKernel_A_gets_A_times_B(bs, pc, pv, multiplier); /* pc= multiplier = pc*inv(diag[row]) */
16419566063dSJacob Faibussowitsch       PetscCall(MatBlockAbs_private(1, bs2, pc, vtmp_abs));
16422efa7f71SHong Zhang       if (vtmp_abs[0] > dt) {         /* apply tolerance dropping rule */
16432efa7f71SHong Zhang         pj = bj + bdiag[row + 1] + 1; /* point to 1st entry of U(row,:) */
16442efa7f71SHong Zhang         pv = ba + bs2 * (bdiag[row + 1] + 1);
16452efa7f71SHong Zhang         nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries in U(row,:), excluding diagonal */
16469371c9d4SSatish Balay         for (j = 0; j < nz; j++) { PetscKernel_A_gets_A_minus_B_times_C(bs, rtmp + bs2 * pj[j], pc, pv + bs2 * j); }
16479566063dSJacob Faibussowitsch         /* PetscCall(PetscLogFlops(bslog*(nz+1.0)-bs)); */
16482efa7f71SHong Zhang       }
16492efa7f71SHong Zhang       row = *bjtmp++;
16502efa7f71SHong Zhang     }
16512efa7f71SHong Zhang 
16522efa7f71SHong Zhang     /* copy sparse rtmp into contiguous vtmp; separate L and U part */
16539371c9d4SSatish Balay     nzi_bl = 0;
16549371c9d4SSatish Balay     j      = 0;
16552efa7f71SHong Zhang     while (jtmp[j] < i) { /* L-part. Note: jtmp is sorted */
16569566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(vtmp + bs2 * j, rtmp + bs2 * jtmp[j], bs2));
16579371c9d4SSatish Balay       nzi_bl++;
16589371c9d4SSatish Balay       j++;
16592efa7f71SHong Zhang     }
16602efa7f71SHong Zhang     nzi_bu = nzi - nzi_bl - 1;
16612efa7f71SHong Zhang 
16622efa7f71SHong Zhang     while (j < nzi) { /* U-part */
16639566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(vtmp + bs2 * j, rtmp + bs2 * jtmp[j], bs2));
16642efa7f71SHong Zhang       j++;
16652efa7f71SHong Zhang     }
16662efa7f71SHong Zhang 
16679566063dSJacob Faibussowitsch     PetscCall(MatBlockAbs_private(nzi, bs2, vtmp, vtmp_abs));
16685f8bbccaSHong Zhang 
16692efa7f71SHong Zhang     bjtmp = bj + bi[i];
16702efa7f71SHong Zhang     batmp = ba + bs2 * bi[i];
16712efa7f71SHong Zhang     /* apply level dropping rule to L part */
16722efa7f71SHong Zhang     ncut  = nzi_al + dtcount;
16732efa7f71SHong Zhang     if (ncut < nzi_bl) {
16749566063dSJacob Faibussowitsch       PetscCall(PetscSortSplitReal(ncut, nzi_bl, vtmp_abs, jtmp));
16759566063dSJacob Faibussowitsch       PetscCall(PetscSortIntWithScalarArray(ncut, jtmp, vtmp));
16762efa7f71SHong Zhang     } else {
16772efa7f71SHong Zhang       ncut = nzi_bl;
16782efa7f71SHong Zhang     }
16792efa7f71SHong Zhang     for (j = 0; j < ncut; j++) {
16802efa7f71SHong Zhang       bjtmp[j] = jtmp[j];
16819566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(batmp + bs2 * j, rtmp + bs2 * bjtmp[j], bs2));
16822efa7f71SHong Zhang     }
16832efa7f71SHong Zhang     bi[i + 1] = bi[i] + ncut;
16842efa7f71SHong Zhang     nzi       = ncut + 1;
16852efa7f71SHong Zhang 
16862efa7f71SHong Zhang     /* apply level dropping rule to U part */
16872efa7f71SHong Zhang     ncut = nzi_au + dtcount;
16882efa7f71SHong Zhang     if (ncut < nzi_bu) {
16899566063dSJacob Faibussowitsch       PetscCall(PetscSortSplitReal(ncut, nzi_bu, vtmp_abs + nzi_bl + 1, jtmp + nzi_bl + 1));
16909566063dSJacob Faibussowitsch       PetscCall(PetscSortIntWithScalarArray(ncut, jtmp + nzi_bl + 1, vtmp + nzi_bl + 1));
16912efa7f71SHong Zhang     } else {
16922efa7f71SHong Zhang       ncut = nzi_bu;
16932efa7f71SHong Zhang     }
16942efa7f71SHong Zhang     nzi += ncut;
16952efa7f71SHong Zhang 
16962efa7f71SHong Zhang     /* mark bdiagonal */
16972efa7f71SHong Zhang     bdiag[i + 1]    = bdiag[i] - (ncut + 1);
16986bce7ff8SHong Zhang     bi[2 * mbs - i] = bi[2 * mbs - i + 1] - (ncut + 1);
16996bce7ff8SHong Zhang 
17002efa7f71SHong Zhang     bjtmp = bj + bdiag[i];
17012efa7f71SHong Zhang     batmp = ba + bs2 * bdiag[i];
17029566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(batmp, rtmp + bs2 * i, bs2));
17032efa7f71SHong Zhang     *bjtmp = i;
17045f8bbccaSHong Zhang 
17052efa7f71SHong Zhang     bjtmp = bj + bdiag[i + 1] + 1;
17062efa7f71SHong Zhang     batmp = ba + (bdiag[i + 1] + 1) * bs2;
17072efa7f71SHong Zhang 
17082efa7f71SHong Zhang     for (k = 0; k < ncut; k++) {
17092efa7f71SHong Zhang       bjtmp[k] = jtmp[nzi_bl + 1 + k];
17109566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(batmp + bs2 * k, rtmp + bs2 * bjtmp[k], bs2));
17112efa7f71SHong Zhang     }
17122efa7f71SHong Zhang 
17132efa7f71SHong Zhang     im[i] = nzi; /* used by PetscLLAddSortedLU() */
17142efa7f71SHong Zhang 
1715a5b23f4aSJose E. Roman     /* invert diagonal block for simpler triangular solves - add shift??? */
17162efa7f71SHong Zhang     batmp = ba + bs2 * bdiag[i];
17175f8bbccaSHong Zhang 
17189566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A(bs, batmp, v_pivots, v_work, allowzeropivot, &zeropivotdetected));
17197b6c816cSBarry Smith     if (zeropivotdetected) B->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
17202efa7f71SHong Zhang   } /* for (i=0; i<mbs; i++) */
17219566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v_work, multiplier, v_pivots));
17222efa7f71SHong Zhang 
17236da515a1SHong Zhang   /* printf("end of L %d, beginning of U %d\n",bi[mbs],bdiag[mbs]); */
172494bad497SJacob Faibussowitsch   PetscCheck(bi[mbs] < bdiag[mbs], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "end of L array %" PetscInt_FMT " cannot >= the beginning of U array %" PetscInt_FMT, bi[mbs], bdiag[mbs]);
17252efa7f71SHong Zhang 
17269566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
17279566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
17282efa7f71SHong Zhang 
17299566063dSJacob Faibussowitsch   PetscCall(PetscLLDestroy(lnk, lnkbt));
17302efa7f71SHong Zhang 
17319566063dSJacob Faibussowitsch   PetscCall(PetscFree2(im, jtmp));
17329566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rtmp, vtmp));
17332efa7f71SHong Zhang 
17349566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(bs2 * B->cmap->n));
17352efa7f71SHong Zhang   b->maxnz = b->nz = bi[mbs] + bdiag[0] - bdiag[mbs];
17362efa7f71SHong Zhang 
17379566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isrow, &row_identity));
17389566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isicol, &icol_identity));
17392efa7f71SHong Zhang   if (row_identity && icol_identity) {
17404dd39f65SShri Abhyankar     B->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering;
17412efa7f71SHong Zhang   } else {
17424dd39f65SShri Abhyankar     B->ops->solve = MatSolve_SeqBAIJ_N;
17432efa7f71SHong Zhang   }
17442efa7f71SHong Zhang 
1745f4259b30SLisandro Dalcin   B->ops->solveadd          = NULL;
1746f4259b30SLisandro Dalcin   B->ops->solvetranspose    = NULL;
1747f4259b30SLisandro Dalcin   B->ops->solvetransposeadd = NULL;
1748f4259b30SLisandro Dalcin   B->ops->matsolve          = NULL;
17492efa7f71SHong Zhang   B->assembled              = PETSC_TRUE;
17502efa7f71SHong Zhang   B->preallocated           = PETSC_TRUE;
1751c8342467SHong Zhang   PetscFunctionReturn(0);
1752c8342467SHong Zhang }
1753