xref: /petsc/src/mat/impls/baij/seq/baijfact.c (revision da81f9329be15cc55f054c8a00978087195c9247)
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 
8d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat B, Mat A, const MatFactorInfo *info)
9d71ae5a4SJacob Faibussowitsch {
10b588c5a2SHong Zhang   Mat             C = B;
11b588c5a2SHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
12b588c5a2SHong Zhang   IS              isrow = b->row, isicol = b->icol;
135a586d82SBarry Smith   const PetscInt *r, *ic;
14bbd65245SShri Abhyankar   PetscInt        i, j, k, nz, nzL, row, *pj;
15bbd65245SShri Abhyankar   const PetscInt  n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j, bs2 = a->bs2;
16bbd65245SShri Abhyankar   const PetscInt *ajtmp, *bjtmp, *bdiag = b->diag;
17bbd65245SShri Abhyankar   MatScalar      *rtmp, *pc, *mwork, *pv;
18bbd65245SShri Abhyankar   MatScalar      *aa = a->a, *v;
19bbd65245SShri Abhyankar   PetscInt        flg;
20182b8fbaSHong Zhang   PetscReal       shift = info->shiftamount;
21a455e926SHong Zhang   PetscBool       allowzeropivot, zeropivotdetected;
22b588c5a2SHong Zhang 
23b588c5a2SHong Zhang   PetscFunctionBegin;
249566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
259566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
260164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
27b588c5a2SHong Zhang 
284c000e2eSHong Zhang   /* generate work space needed by the factorization */
299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork));
309566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(rtmp, bs2 * n));
31b588c5a2SHong Zhang 
32b588c5a2SHong Zhang   for (i = 0; i < n; i++) {
33b588c5a2SHong Zhang     /* zero rtmp */
34b588c5a2SHong Zhang     /* L part */
35b588c5a2SHong Zhang     nz    = bi[i + 1] - bi[i];
36b588c5a2SHong Zhang     bjtmp = bj + bi[i];
3748a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
38b588c5a2SHong Zhang 
39b588c5a2SHong Zhang     /* U part */
400c4413a7SShri Abhyankar     nz    = bdiag[i] - bdiag[i + 1];
410c4413a7SShri Abhyankar     bjtmp = bj + bdiag[i + 1] + 1;
4248a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
430c4413a7SShri Abhyankar 
440c4413a7SShri Abhyankar     /* load in initial (unfactored row) */
450c4413a7SShri Abhyankar     nz    = ai[r[i] + 1] - ai[r[i]];
460c4413a7SShri Abhyankar     ajtmp = aj + ai[r[i]];
470c4413a7SShri Abhyankar     v     = aa + bs2 * ai[r[i]];
4848a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ic[ajtmp[j]], v + bs2 * j, bs2));
490c4413a7SShri Abhyankar 
500c4413a7SShri Abhyankar     /* elimination */
510c4413a7SShri Abhyankar     bjtmp = bj + bi[i];
520c4413a7SShri Abhyankar     nzL   = bi[i + 1] - bi[i];
530c4413a7SShri Abhyankar     for (k = 0; k < nzL; k++) {
540c4413a7SShri Abhyankar       row = bjtmp[k];
550c4413a7SShri Abhyankar       pc  = rtmp + bs2 * row;
56c35f09e5SBarry Smith       for (flg = 0, j = 0; j < bs2; j++) {
57c35f09e5SBarry Smith         if (pc[j] != (PetscScalar)0.0) {
58c35f09e5SBarry Smith           flg = 1;
59c35f09e5SBarry Smith           break;
60c35f09e5SBarry Smith         }
61c35f09e5SBarry Smith       }
620c4413a7SShri Abhyankar       if (flg) {
630c4413a7SShri Abhyankar         pv = b->a + bs2 * bdiag[row];
6496b95a6bSBarry Smith         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
659566063dSJacob Faibussowitsch         PetscCall(PetscKernel_A_gets_A_times_B_2(pc, pv, mwork));
660c4413a7SShri Abhyankar 
67a5b23f4aSJose E. Roman         pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */
680c4413a7SShri Abhyankar         pv = b->a + bs2 * (bdiag[row + 1] + 1);
690c4413a7SShri Abhyankar         nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */
700c4413a7SShri Abhyankar         for (j = 0; j < nz; j++) {
7196b95a6bSBarry Smith           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
720c4413a7SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
730c4413a7SShri Abhyankar           v = rtmp + 4 * pj[j];
749566063dSJacob Faibussowitsch           PetscCall(PetscKernel_A_gets_A_minus_B_times_C_2(v, pc, pv));
750c4413a7SShri Abhyankar           pv += 4;
760c4413a7SShri Abhyankar         }
779566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(16.0 * nz + 12)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
780c4413a7SShri Abhyankar       }
790c4413a7SShri Abhyankar     }
800c4413a7SShri Abhyankar 
810c4413a7SShri Abhyankar     /* finished row so stick it into b->a */
820c4413a7SShri Abhyankar     /* L part */
830c4413a7SShri Abhyankar     pv = b->a + bs2 * bi[i];
840c4413a7SShri Abhyankar     pj = b->j + bi[i];
850c4413a7SShri Abhyankar     nz = bi[i + 1] - bi[i];
8648a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
870c4413a7SShri Abhyankar 
88a5b23f4aSJose E. Roman     /* Mark diagonal and invert diagonal for simpler triangular solves */
890c4413a7SShri Abhyankar     pv = b->a + bs2 * bdiag[i];
900c4413a7SShri Abhyankar     pj = b->j + bdiag[i];
919566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2));
929566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_2(pv, shift, allowzeropivot, &zeropivotdetected));
937b6c816cSBarry Smith     if (zeropivotdetected) B->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
940c4413a7SShri Abhyankar 
950c4413a7SShri Abhyankar     /* U part */
960c4413a7SShri Abhyankar     pv = b->a + bs2 * (bdiag[i + 1] + 1);
970c4413a7SShri Abhyankar     pj = b->j + bdiag[i + 1] + 1;
980c4413a7SShri Abhyankar     nz = bdiag[i] - bdiag[i + 1] - 1;
9948a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
1000c4413a7SShri Abhyankar   }
1010c4413a7SShri Abhyankar 
1029566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rtmp, mwork));
1039566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
1049566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
10526fbe8dcSKarl Rupp 
1064dd39f65SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_2;
1074dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2;
1080c4413a7SShri Abhyankar   C->assembled           = PETSC_TRUE;
10926fbe8dcSKarl Rupp 
1109566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 2 * 2 * 2 * n)); /* from inverting diagonal blocks */
1110c4413a7SShri Abhyankar   PetscFunctionReturn(0);
1120c4413a7SShri Abhyankar }
1130c4413a7SShri Abhyankar 
114d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat B, Mat A, const MatFactorInfo *info)
115d71ae5a4SJacob Faibussowitsch {
1164c000e2eSHong Zhang   Mat             C = B;
1174c000e2eSHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
118bbd65245SShri Abhyankar   PetscInt        i, j, k, nz, nzL, row, *pj;
119bbd65245SShri Abhyankar   const PetscInt  n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j, bs2 = a->bs2;
120bbd65245SShri Abhyankar   const PetscInt *ajtmp, *bjtmp, *bdiag = b->diag;
121bbd65245SShri Abhyankar   MatScalar      *rtmp, *pc, *mwork, *pv;
122bbd65245SShri Abhyankar   MatScalar      *aa = a->a, *v;
123bbd65245SShri Abhyankar   PetscInt        flg;
124182b8fbaSHong Zhang   PetscReal       shift = info->shiftamount;
125a455e926SHong Zhang   PetscBool       allowzeropivot, zeropivotdetected;
1264c000e2eSHong Zhang 
1274c000e2eSHong Zhang   PetscFunctionBegin;
1280164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
1290164db54SHong Zhang 
1304c000e2eSHong Zhang   /* generate work space needed by the factorization */
1319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork));
1329566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(rtmp, bs2 * n));
1334c000e2eSHong Zhang 
1344c000e2eSHong Zhang   for (i = 0; i < n; i++) {
1354c000e2eSHong Zhang     /* zero rtmp */
1364c000e2eSHong Zhang     /* L part */
1374c000e2eSHong Zhang     nz    = bi[i + 1] - bi[i];
1384c000e2eSHong Zhang     bjtmp = bj + bi[i];
13948a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
1404c000e2eSHong Zhang 
1414c000e2eSHong Zhang     /* U part */
142b2b2dd24SShri Abhyankar     nz    = bdiag[i] - bdiag[i + 1];
143b2b2dd24SShri Abhyankar     bjtmp = bj + bdiag[i + 1] + 1;
14448a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
145b2b2dd24SShri Abhyankar 
146b2b2dd24SShri Abhyankar     /* load in initial (unfactored row) */
147b2b2dd24SShri Abhyankar     nz    = ai[i + 1] - ai[i];
148b2b2dd24SShri Abhyankar     ajtmp = aj + ai[i];
149b2b2dd24SShri Abhyankar     v     = aa + bs2 * ai[i];
15048a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ajtmp[j], v + bs2 * j, bs2));
151b2b2dd24SShri Abhyankar 
152b2b2dd24SShri Abhyankar     /* elimination */
153b2b2dd24SShri Abhyankar     bjtmp = bj + bi[i];
154b2b2dd24SShri Abhyankar     nzL   = bi[i + 1] - bi[i];
155b2b2dd24SShri Abhyankar     for (k = 0; k < nzL; k++) {
156b2b2dd24SShri Abhyankar       row = bjtmp[k];
157b2b2dd24SShri Abhyankar       pc  = rtmp + bs2 * row;
158c35f09e5SBarry Smith       for (flg = 0, j = 0; j < bs2; j++) {
159c35f09e5SBarry Smith         if (pc[j] != (PetscScalar)0.0) {
160c35f09e5SBarry Smith           flg = 1;
161c35f09e5SBarry Smith           break;
162c35f09e5SBarry Smith         }
163c35f09e5SBarry Smith       }
164b2b2dd24SShri Abhyankar       if (flg) {
165b2b2dd24SShri Abhyankar         pv = b->a + bs2 * bdiag[row];
16696b95a6bSBarry Smith         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
1679566063dSJacob Faibussowitsch         PetscCall(PetscKernel_A_gets_A_times_B_2(pc, pv, mwork));
168b2b2dd24SShri Abhyankar 
169b2b2dd24SShri Abhyankar         pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */
170b2b2dd24SShri Abhyankar         pv = b->a + bs2 * (bdiag[row + 1] + 1);
171b2b2dd24SShri Abhyankar         nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries in U(row,:) excluding diag */
172b2b2dd24SShri Abhyankar         for (j = 0; j < nz; j++) {
17396b95a6bSBarry Smith           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
174b2b2dd24SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
175b2b2dd24SShri Abhyankar           v = rtmp + 4 * pj[j];
1769566063dSJacob Faibussowitsch           PetscCall(PetscKernel_A_gets_A_minus_B_times_C_2(v, pc, pv));
177b2b2dd24SShri Abhyankar           pv += 4;
178b2b2dd24SShri Abhyankar         }
1799566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(16.0 * nz + 12)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
180b2b2dd24SShri Abhyankar       }
181b2b2dd24SShri Abhyankar     }
182b2b2dd24SShri Abhyankar 
183b2b2dd24SShri Abhyankar     /* finished row so stick it into b->a */
184b2b2dd24SShri Abhyankar     /* L part */
185b2b2dd24SShri Abhyankar     pv = b->a + bs2 * bi[i];
186b2b2dd24SShri Abhyankar     pj = b->j + bi[i];
187b2b2dd24SShri Abhyankar     nz = bi[i + 1] - bi[i];
18848a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
189b2b2dd24SShri Abhyankar 
190a5b23f4aSJose E. Roman     /* Mark diagonal and invert diagonal for simpler triangular solves */
191b2b2dd24SShri Abhyankar     pv = b->a + bs2 * bdiag[i];
192b2b2dd24SShri Abhyankar     pj = b->j + bdiag[i];
1939566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2));
1949566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_2(pv, shift, allowzeropivot, &zeropivotdetected));
1957b6c816cSBarry Smith     if (zeropivotdetected) B->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
196b2b2dd24SShri Abhyankar 
197b2b2dd24SShri Abhyankar     /* U part */
198b2b2dd24SShri Abhyankar     /*
199b2b2dd24SShri Abhyankar     pv = b->a + bs2*bi[2*n-i];
200b2b2dd24SShri Abhyankar     pj = b->j + bi[2*n-i];
201b2b2dd24SShri Abhyankar     nz = bi[2*n-i+1] - bi[2*n-i] - 1;
202b2b2dd24SShri Abhyankar     */
203b2b2dd24SShri Abhyankar     pv = b->a + bs2 * (bdiag[i + 1] + 1);
204b2b2dd24SShri Abhyankar     pj = b->j + bdiag[i + 1] + 1;
205b2b2dd24SShri Abhyankar     nz = bdiag[i] - bdiag[i + 1] - 1;
20648a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
207b2b2dd24SShri Abhyankar   }
2089566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rtmp, mwork));
209ae3d28f0SHong Zhang 
2104dd39f65SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_2_NaturalOrdering;
2119f5c0bcdSBarry Smith   C->ops->forwardsolve   = MatForwardSolve_SeqBAIJ_2_NaturalOrdering;
2129f5c0bcdSBarry Smith   C->ops->backwardsolve  = MatBackwardSolve_SeqBAIJ_2_NaturalOrdering;
2134dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering;
214b2b2dd24SShri Abhyankar   C->assembled           = PETSC_TRUE;
21526fbe8dcSKarl Rupp 
2169566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 2 * 2 * 2 * n)); /* from inverting diagonal blocks */
217b2b2dd24SShri Abhyankar   PetscFunctionReturn(0);
218b2b2dd24SShri Abhyankar }
219b2b2dd24SShri Abhyankar 
220d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_inplace(Mat B, Mat A, const MatFactorInfo *info)
221d71ae5a4SJacob Faibussowitsch {
222719d5645SBarry Smith   Mat             C = B;
2234eeb42bcSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
2247cdf2dbbSSatish Balay   IS              isrow = b->row, isicol = b->icol;
2255d0c19d7SBarry Smith   const PetscInt *r, *ic;
2265d0c19d7SBarry Smith   PetscInt        i, j, n = a->mbs, *bi = b->i, *bj = b->j;
227690b6cddSBarry Smith   PetscInt       *ajtmpold, *ajtmp, nz, row;
228690b6cddSBarry Smith   PetscInt       *diag_offset = b->diag, idx, *ai = a->i, *aj = a->j, *pj;
229329f5518SBarry Smith   MatScalar      *pv, *v, *rtmp, m1, m2, m3, m4, *pc, *w, *x, x1, x2, x3, x4;
2302515f8d2SSatish Balay   MatScalar       p1, p2, p3, p4;
2313f1db9ecSBarry Smith   MatScalar      *ba = b->a, *aa = a->a;
232182b8fbaSHong Zhang   PetscReal       shift = info->shiftamount;
233a455e926SHong Zhang   PetscBool       allowzeropivot, zeropivotdetected;
2344eeb42bcSBarry Smith 
2353a40ed3dSBarry Smith   PetscFunctionBegin;
2360164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
2379566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
2389566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
2399566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(4 * (n + 1), &rtmp));
2404eeb42bcSBarry Smith 
2414eeb42bcSBarry Smith   for (i = 0; i < n; i++) {
2424078e994SBarry Smith     nz    = bi[i + 1] - bi[i];
2434078e994SBarry Smith     ajtmp = bj + bi[i];
2444eeb42bcSBarry Smith     for (j = 0; j < nz; j++) {
2459371c9d4SSatish Balay       x    = rtmp + 4 * ajtmp[j];
2469371c9d4SSatish Balay       x[0] = x[1] = x[2] = x[3] = 0.0;
2474eeb42bcSBarry Smith     }
2484eeb42bcSBarry Smith     /* load in initial (unfactored row) */
2494eeb42bcSBarry Smith     idx      = r[i];
2504078e994SBarry Smith     nz       = ai[idx + 1] - ai[idx];
2514078e994SBarry Smith     ajtmpold = aj + ai[idx];
2524078e994SBarry Smith     v        = aa + 4 * ai[idx];
2534eeb42bcSBarry Smith     for (j = 0; j < nz; j++) {
2544eeb42bcSBarry Smith       x    = rtmp + 4 * ic[ajtmpold[j]];
2559371c9d4SSatish Balay       x[0] = v[0];
2569371c9d4SSatish Balay       x[1] = v[1];
2579371c9d4SSatish Balay       x[2] = v[2];
2589371c9d4SSatish Balay       x[3] = v[3];
2594eeb42bcSBarry Smith       v += 4;
2604eeb42bcSBarry Smith     }
2614eeb42bcSBarry Smith     row = *ajtmp++;
2624eeb42bcSBarry Smith     while (row < i) {
2634eeb42bcSBarry Smith       pc = rtmp + 4 * row;
2649371c9d4SSatish Balay       p1 = pc[0];
2659371c9d4SSatish Balay       p2 = pc[1];
2669371c9d4SSatish Balay       p3 = pc[2];
2679371c9d4SSatish Balay       p4 = pc[3];
268d4a378daSJed Brown       if (p1 != (PetscScalar)0.0 || p2 != (PetscScalar)0.0 || p3 != (PetscScalar)0.0 || p4 != (PetscScalar)0.0) {
2694078e994SBarry Smith         pv    = ba + 4 * diag_offset[row];
2704078e994SBarry Smith         pj    = bj + diag_offset[row] + 1;
2719371c9d4SSatish Balay         x1    = pv[0];
2729371c9d4SSatish Balay         x2    = pv[1];
2739371c9d4SSatish Balay         x3    = pv[2];
2749371c9d4SSatish Balay         x4    = pv[3];
2754eeb42bcSBarry Smith         pc[0] = m1 = p1 * x1 + p3 * x2;
2764eeb42bcSBarry Smith         pc[1] = m2 = p2 * x1 + p4 * x2;
2774eeb42bcSBarry Smith         pc[2] = m3 = p1 * x3 + p3 * x4;
2784eeb42bcSBarry Smith         pc[3] = m4 = p2 * x3 + p4 * x4;
2794078e994SBarry Smith         nz         = bi[row + 1] - diag_offset[row] - 1;
2804eeb42bcSBarry Smith         pv += 4;
2814eeb42bcSBarry Smith         for (j = 0; j < nz; j++) {
2829371c9d4SSatish Balay           x1 = pv[0];
2839371c9d4SSatish Balay           x2 = pv[1];
2849371c9d4SSatish Balay           x3 = pv[2];
2859371c9d4SSatish Balay           x4 = pv[3];
2864eeb42bcSBarry Smith           x  = rtmp + 4 * pj[j];
2874eeb42bcSBarry Smith           x[0] -= m1 * x1 + m3 * x2;
2884eeb42bcSBarry Smith           x[1] -= m2 * x1 + m4 * x2;
2894eeb42bcSBarry Smith           x[2] -= m1 * x3 + m3 * x4;
2904eeb42bcSBarry Smith           x[3] -= m2 * x3 + m4 * x4;
2914eeb42bcSBarry Smith           pv += 4;
2924eeb42bcSBarry Smith         }
2939566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(16.0 * nz + 12.0));
2944eeb42bcSBarry Smith       }
2954eeb42bcSBarry Smith       row = *ajtmp++;
2964eeb42bcSBarry Smith     }
2974eeb42bcSBarry Smith     /* finished row so stick it into b->a */
2984078e994SBarry Smith     pv = ba + 4 * bi[i];
2994078e994SBarry Smith     pj = bj + bi[i];
3004078e994SBarry Smith     nz = bi[i + 1] - bi[i];
3014eeb42bcSBarry Smith     for (j = 0; j < nz; j++) {
3024eeb42bcSBarry Smith       x     = rtmp + 4 * pj[j];
3039371c9d4SSatish Balay       pv[0] = x[0];
3049371c9d4SSatish Balay       pv[1] = x[1];
3059371c9d4SSatish Balay       pv[2] = x[2];
3069371c9d4SSatish Balay       pv[3] = x[3];
3074eeb42bcSBarry Smith       pv += 4;
3084eeb42bcSBarry Smith     }
3094eeb42bcSBarry Smith     /* invert diagonal block */
3104078e994SBarry Smith     w = ba + 4 * diag_offset[i];
3119566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_2(w, shift, allowzeropivot, &zeropivotdetected));
3127b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3134eeb42bcSBarry Smith   }
3144eeb42bcSBarry Smith 
3159566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
3169566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
3179566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
31826fbe8dcSKarl Rupp 
31906e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_2_inplace;
32006e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_inplace;
3214eeb42bcSBarry Smith   C->assembled           = PETSC_TRUE;
32226fbe8dcSKarl Rupp 
3239566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 8 * b->mbs)); /* from inverting diagonal blocks */
3243a40ed3dSBarry Smith   PetscFunctionReturn(0);
3254eeb42bcSBarry Smith }
3264cd38560SBarry Smith /*
3274cd38560SBarry Smith       Version for when blocks are 2 by 2 Using natural ordering
3284cd38560SBarry Smith */
329d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace(Mat C, Mat A, const MatFactorInfo *info)
330d71ae5a4SJacob Faibussowitsch {
3314cd38560SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
332690b6cddSBarry Smith   PetscInt     i, j, n = a->mbs, *bi = b->i, *bj = b->j;
333690b6cddSBarry Smith   PetscInt    *ajtmpold, *ajtmp, nz, row;
334690b6cddSBarry Smith   PetscInt    *diag_offset = b->diag, *ai = a->i, *aj = a->j, *pj;
335329f5518SBarry Smith   MatScalar   *pv, *v, *rtmp, *pc, *w, *x;
3362515f8d2SSatish Balay   MatScalar    p1, p2, p3, p4, m1, m2, m3, m4, x1, x2, x3, x4;
3374cd38560SBarry Smith   MatScalar   *ba = b->a, *aa = a->a;
338182b8fbaSHong Zhang   PetscReal    shift = info->shiftamount;
339a455e926SHong Zhang   PetscBool    allowzeropivot, zeropivotdetected;
3404cd38560SBarry Smith 
3414cd38560SBarry Smith   PetscFunctionBegin;
3420164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
3439566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(4 * (n + 1), &rtmp));
3444cd38560SBarry Smith   for (i = 0; i < n; i++) {
3454cd38560SBarry Smith     nz    = bi[i + 1] - bi[i];
3464cd38560SBarry Smith     ajtmp = bj + bi[i];
3474cd38560SBarry Smith     for (j = 0; j < nz; j++) {
3484cd38560SBarry Smith       x    = rtmp + 4 * ajtmp[j];
3494cd38560SBarry Smith       x[0] = x[1] = x[2] = x[3] = 0.0;
3504cd38560SBarry Smith     }
3514cd38560SBarry Smith     /* load in initial (unfactored row) */
3524cd38560SBarry Smith     nz       = ai[i + 1] - ai[i];
3534cd38560SBarry Smith     ajtmpold = aj + ai[i];
3544cd38560SBarry Smith     v        = aa + 4 * ai[i];
3554cd38560SBarry Smith     for (j = 0; j < nz; j++) {
3564cd38560SBarry Smith       x    = rtmp + 4 * ajtmpold[j];
3579371c9d4SSatish Balay       x[0] = v[0];
3589371c9d4SSatish Balay       x[1] = v[1];
3599371c9d4SSatish Balay       x[2] = v[2];
3609371c9d4SSatish Balay       x[3] = v[3];
3614cd38560SBarry Smith       v += 4;
3624cd38560SBarry Smith     }
3634cd38560SBarry Smith     row = *ajtmp++;
3644cd38560SBarry Smith     while (row < i) {
3654cd38560SBarry Smith       pc = rtmp + 4 * row;
3669371c9d4SSatish Balay       p1 = pc[0];
3679371c9d4SSatish Balay       p2 = pc[1];
3689371c9d4SSatish Balay       p3 = pc[2];
3699371c9d4SSatish Balay       p4 = pc[3];
370d4a378daSJed Brown       if (p1 != (PetscScalar)0.0 || p2 != (PetscScalar)0.0 || p3 != (PetscScalar)0.0 || p4 != (PetscScalar)0.0) {
3714cd38560SBarry Smith         pv    = ba + 4 * diag_offset[row];
3724cd38560SBarry Smith         pj    = bj + diag_offset[row] + 1;
3739371c9d4SSatish Balay         x1    = pv[0];
3749371c9d4SSatish Balay         x2    = pv[1];
3759371c9d4SSatish Balay         x3    = pv[2];
3769371c9d4SSatish Balay         x4    = pv[3];
3774cd38560SBarry Smith         pc[0] = m1 = p1 * x1 + p3 * x2;
3784cd38560SBarry Smith         pc[1] = m2 = p2 * x1 + p4 * x2;
3794cd38560SBarry Smith         pc[2] = m3 = p1 * x3 + p3 * x4;
3804cd38560SBarry Smith         pc[3] = m4 = p2 * x3 + p4 * x4;
3814cd38560SBarry Smith         nz         = bi[row + 1] - diag_offset[row] - 1;
3824cd38560SBarry Smith         pv += 4;
3834cd38560SBarry Smith         for (j = 0; j < nz; j++) {
3849371c9d4SSatish Balay           x1 = pv[0];
3859371c9d4SSatish Balay           x2 = pv[1];
3869371c9d4SSatish Balay           x3 = pv[2];
3879371c9d4SSatish Balay           x4 = pv[3];
3884cd38560SBarry Smith           x  = rtmp + 4 * pj[j];
3894cd38560SBarry Smith           x[0] -= m1 * x1 + m3 * x2;
3904cd38560SBarry Smith           x[1] -= m2 * x1 + m4 * x2;
3914cd38560SBarry Smith           x[2] -= m1 * x3 + m3 * x4;
3924cd38560SBarry Smith           x[3] -= m2 * x3 + m4 * x4;
3934cd38560SBarry Smith           pv += 4;
3944cd38560SBarry Smith         }
3959566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(16.0 * nz + 12.0));
3964cd38560SBarry Smith       }
3974cd38560SBarry Smith       row = *ajtmp++;
3984cd38560SBarry Smith     }
3994cd38560SBarry Smith     /* finished row so stick it into b->a */
4004cd38560SBarry Smith     pv = ba + 4 * bi[i];
4014cd38560SBarry Smith     pj = bj + bi[i];
4024cd38560SBarry Smith     nz = bi[i + 1] - bi[i];
4034cd38560SBarry Smith     for (j = 0; j < nz; j++) {
4044cd38560SBarry Smith       x     = rtmp + 4 * pj[j];
4059371c9d4SSatish Balay       pv[0] = x[0];
4069371c9d4SSatish Balay       pv[1] = x[1];
4079371c9d4SSatish Balay       pv[2] = x[2];
4089371c9d4SSatish Balay       pv[3] = x[3];
4092efa7f71SHong Zhang       /*
4102efa7f71SHong Zhang       printf(" col %d:",pj[j]);
4112efa7f71SHong Zhang       PetscInt j1;
4122efa7f71SHong Zhang       for (j1=0; j1<4; j1++) printf(" %g,",*(pv+j1));
4132efa7f71SHong Zhang       printf("\n");
4142efa7f71SHong Zhang       */
4154cd38560SBarry Smith       pv += 4;
4164cd38560SBarry Smith     }
4174cd38560SBarry Smith     /* invert diagonal block */
4184cd38560SBarry Smith     w = ba + 4 * diag_offset[i];
4199566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_2(w, shift, allowzeropivot, &zeropivotdetected));
4207b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
4214cd38560SBarry Smith   }
4224cd38560SBarry Smith 
4239566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
42426fbe8dcSKarl Rupp 
42506e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_2_NaturalOrdering_inplace;
42606e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering_inplace;
4274cd38560SBarry Smith   C->assembled           = PETSC_TRUE;
42826fbe8dcSKarl Rupp 
4299566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 8 * b->mbs)); /* from inverting diagonal blocks */
4304cd38560SBarry Smith   PetscFunctionReturn(0);
4314cd38560SBarry Smith }
4327fc0212eSBarry Smith 
4337fc0212eSBarry Smith /* ----------------------------------------------------------- */
4347fc0212eSBarry Smith /*
4357fc0212eSBarry Smith      Version for when blocks are 1 by 1.
4367fc0212eSBarry Smith */
437d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat B, Mat A, const MatFactorInfo *info)
438d71ae5a4SJacob Faibussowitsch {
439048b5e81SShri Abhyankar   Mat              C = B;
440048b5e81SShri Abhyankar   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
441048b5e81SShri Abhyankar   IS               isrow = b->row, isicol = b->icol;
442048b5e81SShri Abhyankar   const PetscInt  *r, *ic, *ics;
443048b5e81SShri Abhyankar   const PetscInt   n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j, *bdiag = b->diag;
444048b5e81SShri Abhyankar   PetscInt         i, j, k, nz, nzL, row, *pj;
445048b5e81SShri Abhyankar   const PetscInt  *ajtmp, *bjtmp;
446048b5e81SShri Abhyankar   MatScalar       *rtmp, *pc, multiplier, *pv;
447048b5e81SShri Abhyankar   const MatScalar *aa = a->a, *v;
448ace3abfcSBarry Smith   PetscBool        row_identity, col_identity;
449048b5e81SShri Abhyankar   FactorShiftCtx   sctx;
450048b5e81SShri Abhyankar   const PetscInt  *ddiag;
451048b5e81SShri Abhyankar   PetscReal        rs;
452048b5e81SShri Abhyankar   MatScalar        d;
453048b5e81SShri Abhyankar 
454048b5e81SShri Abhyankar   PetscFunctionBegin;
455048b5e81SShri Abhyankar   /* MatPivotSetUp(): initialize shift context sctx */
4569566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&sctx, sizeof(FactorShiftCtx)));
457048b5e81SShri Abhyankar 
458048b5e81SShri Abhyankar   if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
459048b5e81SShri Abhyankar     ddiag          = a->diag;
460048b5e81SShri Abhyankar     sctx.shift_top = info->zeropivot;
461048b5e81SShri Abhyankar     for (i = 0; i < n; i++) {
462048b5e81SShri Abhyankar       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
463048b5e81SShri Abhyankar       d  = (aa)[ddiag[i]];
464048b5e81SShri Abhyankar       rs = -PetscAbsScalar(d) - PetscRealPart(d);
465048b5e81SShri Abhyankar       v  = aa + ai[i];
466048b5e81SShri Abhyankar       nz = ai[i + 1] - ai[i];
46726fbe8dcSKarl Rupp       for (j = 0; j < nz; j++) rs += PetscAbsScalar(v[j]);
468048b5e81SShri Abhyankar       if (rs > sctx.shift_top) sctx.shift_top = rs;
469048b5e81SShri Abhyankar     }
470048b5e81SShri Abhyankar     sctx.shift_top *= 1.1;
471048b5e81SShri Abhyankar     sctx.nshift_max = 5;
472048b5e81SShri Abhyankar     sctx.shift_lo   = 0.;
473048b5e81SShri Abhyankar     sctx.shift_hi   = 1.;
474048b5e81SShri Abhyankar   }
475048b5e81SShri Abhyankar 
4769566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
4779566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
4789566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n + 1, &rtmp));
479048b5e81SShri Abhyankar   ics = ic;
480048b5e81SShri Abhyankar 
481048b5e81SShri Abhyankar   do {
482048b5e81SShri Abhyankar     sctx.newshift = PETSC_FALSE;
483048b5e81SShri Abhyankar     for (i = 0; i < n; i++) {
484048b5e81SShri Abhyankar       /* zero rtmp */
485048b5e81SShri Abhyankar       /* L part */
486048b5e81SShri Abhyankar       nz    = bi[i + 1] - bi[i];
487048b5e81SShri Abhyankar       bjtmp = bj + bi[i];
488048b5e81SShri Abhyankar       for (j = 0; j < nz; j++) rtmp[bjtmp[j]] = 0.0;
489048b5e81SShri Abhyankar 
490048b5e81SShri Abhyankar       /* U part */
491048b5e81SShri Abhyankar       nz    = bdiag[i] - bdiag[i + 1];
492048b5e81SShri Abhyankar       bjtmp = bj + bdiag[i + 1] + 1;
493048b5e81SShri Abhyankar       for (j = 0; j < nz; j++) rtmp[bjtmp[j]] = 0.0;
494048b5e81SShri Abhyankar 
495048b5e81SShri Abhyankar       /* load in initial (unfactored row) */
496048b5e81SShri Abhyankar       nz    = ai[r[i] + 1] - ai[r[i]];
497048b5e81SShri Abhyankar       ajtmp = aj + ai[r[i]];
498048b5e81SShri Abhyankar       v     = aa + ai[r[i]];
49926fbe8dcSKarl Rupp       for (j = 0; j < nz; j++) rtmp[ics[ajtmp[j]]] = v[j];
50026fbe8dcSKarl Rupp 
501048b5e81SShri Abhyankar       /* ZeropivotApply() */
502048b5e81SShri Abhyankar       rtmp[i] += sctx.shift_amount; /* shift the diagonal of the matrix */
503048b5e81SShri Abhyankar 
504048b5e81SShri Abhyankar       /* elimination */
505048b5e81SShri Abhyankar       bjtmp = bj + bi[i];
506048b5e81SShri Abhyankar       row   = *bjtmp++;
507048b5e81SShri Abhyankar       nzL   = bi[i + 1] - bi[i];
508048b5e81SShri Abhyankar       for (k = 0; k < nzL; k++) {
509048b5e81SShri Abhyankar         pc = rtmp + row;
510d4a378daSJed Brown         if (*pc != (PetscScalar)0.0) {
511048b5e81SShri Abhyankar           pv         = b->a + bdiag[row];
512048b5e81SShri Abhyankar           multiplier = *pc * (*pv);
513048b5e81SShri Abhyankar           *pc        = multiplier;
51426fbe8dcSKarl Rupp 
515048b5e81SShri Abhyankar           pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */
516048b5e81SShri Abhyankar           pv = b->a + bdiag[row + 1] + 1;
517048b5e81SShri Abhyankar           nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries in U(row,:) excluding diag */
518048b5e81SShri Abhyankar           for (j = 0; j < nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
5199566063dSJacob Faibussowitsch           PetscCall(PetscLogFlops(2.0 * nz));
520048b5e81SShri Abhyankar         }
521048b5e81SShri Abhyankar         row = *bjtmp++;
522048b5e81SShri Abhyankar       }
523048b5e81SShri Abhyankar 
524048b5e81SShri Abhyankar       /* finished row so stick it into b->a */
525048b5e81SShri Abhyankar       rs = 0.0;
526048b5e81SShri Abhyankar       /* L part */
527048b5e81SShri Abhyankar       pv = b->a + bi[i];
528048b5e81SShri Abhyankar       pj = b->j + bi[i];
529048b5e81SShri Abhyankar       nz = bi[i + 1] - bi[i];
530048b5e81SShri Abhyankar       for (j = 0; j < nz; j++) {
5319371c9d4SSatish Balay         pv[j] = rtmp[pj[j]];
5329371c9d4SSatish Balay         rs += PetscAbsScalar(pv[j]);
533048b5e81SShri Abhyankar       }
534048b5e81SShri Abhyankar 
535048b5e81SShri Abhyankar       /* U part */
536048b5e81SShri Abhyankar       pv = b->a + bdiag[i + 1] + 1;
537048b5e81SShri Abhyankar       pj = b->j + bdiag[i + 1] + 1;
538048b5e81SShri Abhyankar       nz = bdiag[i] - bdiag[i + 1] - 1;
539048b5e81SShri Abhyankar       for (j = 0; j < nz; j++) {
5409371c9d4SSatish Balay         pv[j] = rtmp[pj[j]];
5419371c9d4SSatish Balay         rs += PetscAbsScalar(pv[j]);
542048b5e81SShri Abhyankar       }
543048b5e81SShri Abhyankar 
544048b5e81SShri Abhyankar       sctx.rs = rs;
545048b5e81SShri Abhyankar       sctx.pv = rtmp[i];
5469566063dSJacob Faibussowitsch       PetscCall(MatPivotCheck(B, A, info, &sctx, i));
547048b5e81SShri Abhyankar       if (sctx.newshift) break; /* break for-loop */
548048b5e81SShri Abhyankar       rtmp[i] = sctx.pv;        /* sctx.pv might be updated in the case of MAT_SHIFT_INBLOCKS */
549048b5e81SShri Abhyankar 
550a5b23f4aSJose E. Roman       /* Mark diagonal and invert diagonal for simpler triangular solves */
551048b5e81SShri Abhyankar       pv  = b->a + bdiag[i];
552d4a378daSJed Brown       *pv = (PetscScalar)1.0 / rtmp[i];
553048b5e81SShri Abhyankar 
554048b5e81SShri Abhyankar     } /* endof for (i=0; i<n; i++) { */
555048b5e81SShri Abhyankar 
556048b5e81SShri Abhyankar     /* MatPivotRefine() */
557048b5e81SShri Abhyankar     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE && !sctx.newshift && sctx.shift_fraction > 0 && sctx.nshift < sctx.nshift_max) {
558048b5e81SShri Abhyankar       /*
559048b5e81SShri Abhyankar        * if no shift in this attempt & shifting & started shifting & can refine,
560048b5e81SShri Abhyankar        * then try lower shift
561048b5e81SShri Abhyankar        */
562048b5e81SShri Abhyankar       sctx.shift_hi       = sctx.shift_fraction;
563048b5e81SShri Abhyankar       sctx.shift_fraction = (sctx.shift_hi + sctx.shift_lo) / 2.;
564048b5e81SShri Abhyankar       sctx.shift_amount   = sctx.shift_fraction * sctx.shift_top;
565048b5e81SShri Abhyankar       sctx.newshift       = PETSC_TRUE;
566048b5e81SShri Abhyankar       sctx.nshift++;
567048b5e81SShri Abhyankar     }
568048b5e81SShri Abhyankar   } while (sctx.newshift);
569048b5e81SShri Abhyankar 
5709566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
5719566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
5729566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
573048b5e81SShri Abhyankar 
5749566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isrow, &row_identity));
5759566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isicol, &col_identity));
576048b5e81SShri Abhyankar   if (row_identity && col_identity) {
577048b5e81SShri Abhyankar     C->ops->solve          = MatSolve_SeqBAIJ_1_NaturalOrdering;
5789f5c0bcdSBarry Smith     C->ops->forwardsolve   = MatForwardSolve_SeqBAIJ_1_NaturalOrdering;
5799f5c0bcdSBarry Smith     C->ops->backwardsolve  = MatBackwardSolve_SeqBAIJ_1_NaturalOrdering;
58093fd935bSShri Abhyankar     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering;
581048b5e81SShri Abhyankar   } else {
582048b5e81SShri Abhyankar     C->ops->solve          = MatSolve_SeqBAIJ_1;
58393fd935bSShri Abhyankar     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1;
584048b5e81SShri Abhyankar   }
585048b5e81SShri Abhyankar   C->assembled = PETSC_TRUE;
5869566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(C->cmap->n));
587048b5e81SShri Abhyankar 
588048b5e81SShri Abhyankar   /* MatShiftView(A,info,&sctx) */
589048b5e81SShri Abhyankar   if (sctx.nshift) {
590048b5e81SShri Abhyankar     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
5919566063dSJacob 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));
592048b5e81SShri Abhyankar     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
5939566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shift_nz tries %" PetscInt_FMT ", shift_amount %g\n", sctx.nshift, (double)sctx.shift_amount));
594048b5e81SShri Abhyankar     } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) {
5959566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shift_inblocks applied %" PetscInt_FMT ", each shift_amount %g\n", sctx.nshift, (double)info->shiftamount));
596048b5e81SShri Abhyankar     }
597048b5e81SShri Abhyankar   }
598048b5e81SShri Abhyankar   PetscFunctionReturn(0);
599048b5e81SShri Abhyankar }
600048b5e81SShri Abhyankar 
601d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1_inplace(Mat C, Mat A, const MatFactorInfo *info)
602d71ae5a4SJacob Faibussowitsch {
6037fc0212eSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
6047cdf2dbbSSatish Balay   IS              isrow = b->row, isicol = b->icol;
6055d0c19d7SBarry Smith   const PetscInt *r, *ic;
6065d0c19d7SBarry Smith   PetscInt        i, j, n = a->mbs, *bi = b->i, *bj = b->j;
607690b6cddSBarry Smith   PetscInt       *ajtmpold, *ajtmp, nz, row, *ai = a->i, *aj = a->j;
608690b6cddSBarry Smith   PetscInt       *diag_offset = b->diag, diag, *pj;
609329f5518SBarry Smith   MatScalar      *pv, *v, *rtmp, multiplier, *pc;
6103f1db9ecSBarry Smith   MatScalar      *ba = b->a, *aa = a->a;
611ace3abfcSBarry Smith   PetscBool       row_identity, col_identity;
6127fc0212eSBarry Smith 
6133a40ed3dSBarry Smith   PetscFunctionBegin;
6149566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
6159566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
6169566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n + 1, &rtmp));
6177fc0212eSBarry Smith 
6187fc0212eSBarry Smith   for (i = 0; i < n; i++) {
6194078e994SBarry Smith     nz    = bi[i + 1] - bi[i];
6204078e994SBarry Smith     ajtmp = bj + bi[i];
6217fc0212eSBarry Smith     for (j = 0; j < nz; j++) rtmp[ajtmp[j]] = 0.0;
6227fc0212eSBarry Smith 
6237fc0212eSBarry Smith     /* load in initial (unfactored row) */
6244078e994SBarry Smith     nz       = ai[r[i] + 1] - ai[r[i]];
6254078e994SBarry Smith     ajtmpold = aj + ai[r[i]];
6264078e994SBarry Smith     v        = aa + ai[r[i]];
6277fc0212eSBarry Smith     for (j = 0; j < nz; j++) rtmp[ic[ajtmpold[j]]] = v[j];
6287fc0212eSBarry Smith 
6297fc0212eSBarry Smith     row = *ajtmp++;
6307fc0212eSBarry Smith     while (row < i) {
6317fc0212eSBarry Smith       pc = rtmp + row;
6327fc0212eSBarry Smith       if (*pc != 0.0) {
6334078e994SBarry Smith         pv         = ba + diag_offset[row];
6344078e994SBarry Smith         pj         = bj + diag_offset[row] + 1;
6357fc0212eSBarry Smith         multiplier = *pc * *pv++;
6367fc0212eSBarry Smith         *pc        = multiplier;
6374078e994SBarry Smith         nz         = bi[row + 1] - diag_offset[row] - 1;
6387fc0212eSBarry Smith         for (j = 0; j < nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
6399566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(1.0 + 2.0 * nz));
6407fc0212eSBarry Smith       }
6417fc0212eSBarry Smith       row = *ajtmp++;
6427fc0212eSBarry Smith     }
6437fc0212eSBarry Smith     /* finished row so stick it into b->a */
6444078e994SBarry Smith     pv = ba + bi[i];
6454078e994SBarry Smith     pj = bj + bi[i];
6464078e994SBarry Smith     nz = bi[i + 1] - bi[i];
64726fbe8dcSKarl Rupp     for (j = 0; j < nz; j++) pv[j] = rtmp[pj[j]];
6484078e994SBarry Smith     diag = diag_offset[i] - bi[i];
6497fc0212eSBarry Smith     /* check pivot entry for current row */
65008401ef6SPierre 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);
6517fc0212eSBarry Smith     pv[diag] = 1.0 / pv[diag];
6527fc0212eSBarry Smith   }
6537fc0212eSBarry Smith 
6549566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
6559566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
6569566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
6579566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isrow, &row_identity));
6589566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isicol, &col_identity));
659f58c8c74SBarry Smith   if (row_identity && col_identity) {
66006e38f1dSHong Zhang     C->ops->solve          = MatSolve_SeqBAIJ_1_NaturalOrdering_inplace;
66106e38f1dSHong Zhang     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering_inplace;
662f58c8c74SBarry Smith   } else {
66306e38f1dSHong Zhang     C->ops->solve          = MatSolve_SeqBAIJ_1_inplace;
66406e38f1dSHong Zhang     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_inplace;
665f58c8c74SBarry Smith   }
6667fc0212eSBarry Smith   C->assembled = PETSC_TRUE;
6679566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(C->cmap->n));
6683a40ed3dSBarry Smith   PetscFunctionReturn(0);
6697fc0212eSBarry Smith }
6707fc0212eSBarry Smith 
671d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_petsc(Mat A, MatSolverType *type)
672d71ae5a4SJacob Faibussowitsch {
6734ac6704cSBarry Smith   PetscFunctionBegin;
6744ac6704cSBarry Smith   *type = MATSOLVERPETSC;
6754ac6704cSBarry Smith   PetscFunctionReturn(0);
6764ac6704cSBarry Smith }
6774ac6704cSBarry Smith 
678d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A, MatFactorType ftype, Mat *B)
679d71ae5a4SJacob Faibussowitsch {
680d0f46423SBarry Smith   PetscInt n = A->rmap->n;
681b24902e0SBarry Smith 
682b24902e0SBarry Smith   PetscFunctionBegin;
683b499a2aeSHong Zhang #if defined(PETSC_USE_COMPLEX)
684b94d7dedSBarry 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");
685b499a2aeSHong Zhang #endif
6869566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
6879566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*B, n, n, n, n));
688c8342467SHong Zhang   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) {
6899566063dSJacob Faibussowitsch     PetscCall(MatSetType(*B, MATSEQBAIJ));
69026fbe8dcSKarl Rupp 
6914dd39f65SShri Abhyankar     (*B)->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqBAIJ;
6924dd39f65SShri Abhyankar     (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqBAIJ;
6939566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGND, (char **)&(*B)->preferredordering[MAT_FACTOR_LU]));
6949566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGNATURAL, (char **)&(*B)->preferredordering[MAT_FACTOR_ILU]));
6959566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGNATURAL, (char **)&(*B)->preferredordering[MAT_FACTOR_ILUDT]));
696b24902e0SBarry Smith   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
6979566063dSJacob Faibussowitsch     PetscCall(MatSetType(*B, MATSEQSBAIJ));
6989566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(*B, A->rmap->bs, MAT_SKIP_ALLOCATION, NULL));
69926fbe8dcSKarl Rupp 
7005c9eb25fSBarry Smith     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqBAIJ;
7015c9eb25fSBarry Smith     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ;
7024ac6704cSBarry Smith     /*  Future optimization would be direct symbolic and numerical factorization for BAIJ to support orderings and Cholesky, instead of first converting to SBAIJ */
7039566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGNATURAL, (char **)&(*B)->preferredordering[MAT_FACTOR_CHOLESKY]));
7049566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGNATURAL, (char **)&(*B)->preferredordering[MAT_FACTOR_ICC]));
705e32f2f54SBarry Smith   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported");
706d5f3da31SBarry Smith   (*B)->factortype     = ftype;
707f73b0415SBarry Smith   (*B)->canuseordering = PETSC_TRUE;
70800c67f3bSHong Zhang 
7099566063dSJacob Faibussowitsch   PetscCall(PetscFree((*B)->solvertype));
7109566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERPETSC, &(*B)->solvertype));
7119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)*B, "MatFactorGetSolverType_C", MatFactorGetSolverType_petsc));
712b24902e0SBarry Smith   PetscFunctionReturn(0);
713b24902e0SBarry Smith }
714273d9f13SBarry Smith 
7155d517e7eSBarry Smith /* ----------------------------------------------------------- */
716d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactor_SeqBAIJ(Mat A, IS row, IS col, const MatFactorInfo *info)
717d71ae5a4SJacob Faibussowitsch {
7185d517e7eSBarry Smith   Mat C;
7195d517e7eSBarry Smith 
7203a40ed3dSBarry Smith   PetscFunctionBegin;
7219566063dSJacob Faibussowitsch   PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_LU, &C));
7229566063dSJacob Faibussowitsch   PetscCall(MatLUFactorSymbolic(C, A, row, col, info));
7239566063dSJacob Faibussowitsch   PetscCall(MatLUFactorNumeric(C, A, info));
72426fbe8dcSKarl Rupp 
725db4efbfdSBarry Smith   A->ops->solve          = C->ops->solve;
726db4efbfdSBarry Smith   A->ops->solvetranspose = C->ops->solvetranspose;
72726fbe8dcSKarl Rupp 
7289566063dSJacob Faibussowitsch   PetscCall(MatHeaderMerge(A, &C));
7293a40ed3dSBarry Smith   PetscFunctionReturn(0);
7305d517e7eSBarry Smith }
7314cd38560SBarry Smith 
732c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h>
733d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat C, Mat A, const MatFactorInfo *info)
734d71ae5a4SJacob Faibussowitsch {
73578910aadSHong Zhang   Mat_SeqBAIJ    *a  = (Mat_SeqBAIJ *)A->data;
73678910aadSHong Zhang   Mat_SeqSBAIJ   *b  = (Mat_SeqSBAIJ *)C->data;
73778910aadSHong Zhang   IS              ip = b->row;
7385d0c19d7SBarry Smith   const PetscInt *rip;
7395d0c19d7SBarry Smith   PetscInt        i, j, mbs = a->mbs, bs = A->rmap->bs, *bi = b->i, *bj = b->j, *bcol;
74078910aadSHong Zhang   PetscInt       *ai = a->i, *aj = a->j;
74178910aadSHong Zhang   PetscInt        k, jmin, jmax, *jl, *il, col, nexti, ili, nz;
74278910aadSHong Zhang   MatScalar      *rtmp, *ba = b->a, *bval, *aa = a->a, dk, uikdi;
74307b50cabSHong Zhang   PetscReal       rs;
7440e95ead3SHong Zhang   FactorShiftCtx  sctx;
74578910aadSHong Zhang 
746c05c3958SHong Zhang   PetscFunctionBegin;
74707b50cabSHong Zhang   if (bs > 1) { /* convert A to a SBAIJ matrix and apply Cholesky factorization from it */
74848a46eb9SPierre Jolivet     if (!a->sbaijMat) PetscCall(MatConvert(A, MATSEQSBAIJ, MAT_INITIAL_MATRIX, &a->sbaijMat));
7499566063dSJacob Faibussowitsch     PetscCall((a->sbaijMat)->ops->choleskyfactornumeric(C, a->sbaijMat, info));
7509566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&a->sbaijMat));
7516ad2eaddSHong Zhang     PetscFunctionReturn(0);
7526ad2eaddSHong Zhang   }
75378910aadSHong Zhang 
75407b50cabSHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
7559566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&sctx, sizeof(FactorShiftCtx)));
7563cea4cbeSHong Zhang 
7579566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(ip, &rip));
7589566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(mbs, &rtmp, mbs, &il, mbs, &jl));
75978910aadSHong Zhang 
76075567043SBarry Smith   sctx.shift_amount = 0.;
7613cea4cbeSHong Zhang   sctx.nshift       = 0;
76278910aadSHong Zhang   do {
76307b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
76478910aadSHong Zhang     for (i = 0; i < mbs; i++) {
7659371c9d4SSatish Balay       rtmp[i] = 0.0;
7669371c9d4SSatish Balay       jl[i]   = mbs;
7679371c9d4SSatish Balay       il[0]   = 0;
76878910aadSHong Zhang     }
76978910aadSHong Zhang 
77078910aadSHong Zhang     for (k = 0; k < mbs; k++) {
77178910aadSHong Zhang       bval = ba + bi[k];
77278910aadSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
7739371c9d4SSatish Balay       jmin = ai[rip[k]];
7749371c9d4SSatish Balay       jmax = ai[rip[k] + 1];
77578910aadSHong Zhang       for (j = jmin; j < jmax; j++) {
77678910aadSHong Zhang         col = rip[aj[j]];
77778910aadSHong Zhang         if (col >= k) { /* only take upper triangular entry */
77878910aadSHong Zhang           rtmp[col] = aa[j];
77978910aadSHong Zhang           *bval++   = 0.0; /* for in-place factorization */
78078910aadSHong Zhang         }
78178910aadSHong Zhang       }
7823cea4cbeSHong Zhang 
7833cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
7843cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
78578910aadSHong Zhang 
78678910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
78778910aadSHong Zhang       dk = rtmp[k];
78878910aadSHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
78978910aadSHong Zhang 
79078910aadSHong Zhang       while (i < k) {
79178910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
79278910aadSHong Zhang 
79378910aadSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
79478910aadSHong Zhang         ili   = il[i];                /* index of first nonzero element in U(i,k:bms-1) */
79578910aadSHong Zhang         uikdi = -ba[ili] * ba[bi[i]]; /* diagonal(k) */
79678910aadSHong Zhang         dk += uikdi * ba[ili];
79778910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
79878910aadSHong Zhang 
79978910aadSHong Zhang         /* add multiple of row i to k-th row */
8009371c9d4SSatish Balay         jmin = ili + 1;
8019371c9d4SSatish Balay         jmax = bi[i + 1];
80278910aadSHong Zhang         if (jmin < jmax) {
80378910aadSHong Zhang           for (j = jmin; j < jmax; j++) rtmp[bj[j]] += uikdi * ba[j];
80478910aadSHong Zhang           /* update il and jl for row i */
80578910aadSHong Zhang           il[i] = jmin;
8069371c9d4SSatish Balay           j     = bj[jmin];
8079371c9d4SSatish Balay           jl[i] = jl[j];
8089371c9d4SSatish Balay           jl[j] = i;
80978910aadSHong Zhang         }
81078910aadSHong Zhang         i = nexti;
81178910aadSHong Zhang       }
81278910aadSHong Zhang 
8133cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
8143cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
8153cea4cbeSHong Zhang       rs   = 0.0;
81678910aadSHong Zhang       jmin = bi[k] + 1;
81778910aadSHong Zhang       nz   = bi[k + 1] - jmin;
81878910aadSHong Zhang       if (nz) {
81978910aadSHong Zhang         bcol = bj + jmin;
82078910aadSHong Zhang         while (nz--) {
82189f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
82289f845c8SHong Zhang           bcol++;
82378910aadSHong Zhang         }
82478910aadSHong Zhang       }
8253cea4cbeSHong Zhang 
8263cea4cbeSHong Zhang       sctx.rs = rs;
8273cea4cbeSHong Zhang       sctx.pv = dk;
8289566063dSJacob Faibussowitsch       PetscCall(MatPivotCheck(C, A, info, &sctx, k));
82907b50cabSHong Zhang       if (sctx.newshift) break;
8300e95ead3SHong Zhang       dk = sctx.pv;
83178910aadSHong Zhang 
83278910aadSHong Zhang       /* copy data into U(k,:) */
83378910aadSHong Zhang       ba[bi[k]] = 1.0 / dk; /* U(k,k) */
8349371c9d4SSatish Balay       jmin      = bi[k] + 1;
8359371c9d4SSatish Balay       jmax      = bi[k + 1];
83678910aadSHong Zhang       if (jmin < jmax) {
83778910aadSHong Zhang         for (j = jmin; j < jmax; j++) {
8389371c9d4SSatish Balay           col       = bj[j];
8399371c9d4SSatish Balay           ba[j]     = rtmp[col];
8409371c9d4SSatish Balay           rtmp[col] = 0.0;
84178910aadSHong Zhang         }
84278910aadSHong Zhang         /* add the k-th row into il and jl */
84378910aadSHong Zhang         il[k] = jmin;
8449371c9d4SSatish Balay         i     = bj[jmin];
8459371c9d4SSatish Balay         jl[k] = jl[i];
8469371c9d4SSatish Balay         jl[i] = k;
84778910aadSHong Zhang       }
84878910aadSHong Zhang     }
84907b50cabSHong Zhang   } while (sctx.newshift);
8509566063dSJacob Faibussowitsch   PetscCall(PetscFree3(rtmp, il, jl));
85178910aadSHong Zhang 
8529566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(ip, &rip));
85326fbe8dcSKarl Rupp 
85478910aadSHong Zhang   C->assembled    = PETSC_TRUE;
85578910aadSHong Zhang   C->preallocated = PETSC_TRUE;
85626fbe8dcSKarl Rupp 
8579566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(C->rmap->N));
8583cea4cbeSHong Zhang   if (sctx.nshift) {
859f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
8609566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shiftpd tries %" PetscInt_FMT ", shift_amount %g\n", sctx.nshift, (double)sctx.shift_amount));
861f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
8629566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shiftnz tries %" PetscInt_FMT ", shift_amount %g\n", sctx.nshift, (double)sctx.shift_amount));
86378910aadSHong Zhang     }
86478910aadSHong Zhang   }
865c05c3958SHong Zhang   PetscFunctionReturn(0);
866c05c3958SHong Zhang }
8674cd38560SBarry Smith 
868d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat C, Mat A, const MatFactorInfo *info)
869d71ae5a4SJacob Faibussowitsch {
87078910aadSHong Zhang   Mat_SeqBAIJ   *a = (Mat_SeqBAIJ *)A->data;
87178910aadSHong Zhang   Mat_SeqSBAIJ  *b = (Mat_SeqSBAIJ *)C->data;
8723cea4cbeSHong Zhang   PetscInt       i, j, am = a->mbs;
87378910aadSHong Zhang   PetscInt      *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j;
8743cea4cbeSHong Zhang   PetscInt       k, jmin, *jl, *il, nexti, ili, *acol, *bcol, nz;
87578910aadSHong Zhang   MatScalar     *rtmp, *ba = b->a, *aa = a->a, dk, uikdi, *aval, *bval;
8760e95ead3SHong Zhang   PetscReal      rs;
8770e95ead3SHong Zhang   FactorShiftCtx sctx;
87878910aadSHong Zhang 
879c05c3958SHong Zhang   PetscFunctionBegin;
8800e95ead3SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
8819566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&sctx, sizeof(FactorShiftCtx)));
8823cea4cbeSHong Zhang 
8839566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(am, &rtmp, am, &il, am, &jl));
88478910aadSHong Zhang 
88578910aadSHong Zhang   do {
88607b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
88778910aadSHong Zhang     for (i = 0; i < am; i++) {
8889371c9d4SSatish Balay       rtmp[i] = 0.0;
8899371c9d4SSatish Balay       jl[i]   = am;
8909371c9d4SSatish Balay       il[0]   = 0;
89178910aadSHong Zhang     }
89278910aadSHong Zhang 
89378910aadSHong Zhang     for (k = 0; k < am; k++) {
89478910aadSHong Zhang       /* initialize k-th row with elements nonzero in row perm(k) of A */
89578910aadSHong Zhang       nz   = ai[k + 1] - ai[k];
89678910aadSHong Zhang       acol = aj + ai[k];
89778910aadSHong Zhang       aval = aa + ai[k];
89878910aadSHong Zhang       bval = ba + bi[k];
89978910aadSHong Zhang       while (nz--) {
90078910aadSHong Zhang         if (*acol < k) { /* skip lower triangular entries */
9019371c9d4SSatish Balay           acol++;
9029371c9d4SSatish Balay           aval++;
90378910aadSHong Zhang         } else {
90478910aadSHong Zhang           rtmp[*acol++] = *aval++;
90578910aadSHong Zhang           *bval++       = 0.0; /* for in-place factorization */
90678910aadSHong Zhang         }
90778910aadSHong Zhang       }
9083cea4cbeSHong Zhang 
9093cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
9103cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
91178910aadSHong Zhang 
91278910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
91378910aadSHong Zhang       dk = rtmp[k];
91478910aadSHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
91578910aadSHong Zhang 
91678910aadSHong Zhang       while (i < k) {
91778910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
91878910aadSHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
91978910aadSHong Zhang         ili   = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
92078910aadSHong Zhang         uikdi = -ba[ili] * ba[bi[i]];
92178910aadSHong Zhang         dk += uikdi * ba[ili];
92278910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
92378910aadSHong Zhang 
92478910aadSHong Zhang         /* add multiple of row i to k-th row ... */
92578910aadSHong Zhang         jmin = ili + 1;
92678910aadSHong Zhang         nz   = bi[i + 1] - jmin;
92778910aadSHong Zhang         if (nz > 0) {
92878910aadSHong Zhang           bcol = bj + jmin;
92978910aadSHong Zhang           bval = ba + jmin;
93078910aadSHong Zhang           while (nz--) rtmp[*bcol++] += uikdi * (*bval++);
93178910aadSHong Zhang           /* update il and jl for i-th row */
93278910aadSHong Zhang           il[i] = jmin;
9339371c9d4SSatish Balay           j     = bj[jmin];
9349371c9d4SSatish Balay           jl[i] = jl[j];
9359371c9d4SSatish Balay           jl[j] = i;
93678910aadSHong Zhang         }
93778910aadSHong Zhang         i = nexti;
93878910aadSHong Zhang       }
93978910aadSHong Zhang 
9403cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
9413cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
9423cea4cbeSHong Zhang       rs   = 0.0;
94378910aadSHong Zhang       jmin = bi[k] + 1;
94478910aadSHong Zhang       nz   = bi[k + 1] - jmin;
94578910aadSHong Zhang       if (nz) {
94678910aadSHong Zhang         bcol = bj + jmin;
94778910aadSHong Zhang         while (nz--) {
94889f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
94989f845c8SHong Zhang           bcol++;
95078910aadSHong Zhang         }
95178910aadSHong Zhang       }
9523cea4cbeSHong Zhang 
9533cea4cbeSHong Zhang       sctx.rs = rs;
9543cea4cbeSHong Zhang       sctx.pv = dk;
9559566063dSJacob Faibussowitsch       PetscCall(MatPivotCheck(C, A, info, &sctx, k));
95607b50cabSHong Zhang       if (sctx.newshift) break; /* sctx.shift_amount is updated */
9570e95ead3SHong Zhang       dk = sctx.pv;
95878910aadSHong Zhang 
95978910aadSHong Zhang       /* copy data into U(k,:) */
96078910aadSHong Zhang       ba[bi[k]] = 1.0 / dk;
96178910aadSHong Zhang       jmin      = bi[k] + 1;
96278910aadSHong Zhang       nz        = bi[k + 1] - jmin;
96378910aadSHong Zhang       if (nz) {
96478910aadSHong Zhang         bcol = bj + jmin;
96578910aadSHong Zhang         bval = ba + jmin;
96678910aadSHong Zhang         while (nz--) {
96778910aadSHong Zhang           *bval++       = rtmp[*bcol];
96878910aadSHong Zhang           rtmp[*bcol++] = 0.0;
96978910aadSHong Zhang         }
97078910aadSHong Zhang         /* add k-th row into il and jl */
97178910aadSHong Zhang         il[k] = jmin;
9729371c9d4SSatish Balay         i     = bj[jmin];
9739371c9d4SSatish Balay         jl[k] = jl[i];
9749371c9d4SSatish Balay         jl[i] = k;
97578910aadSHong Zhang       }
97678910aadSHong Zhang     }
97707b50cabSHong Zhang   } while (sctx.newshift);
9789566063dSJacob Faibussowitsch   PetscCall(PetscFree3(rtmp, il, jl));
97978910aadSHong Zhang 
9800a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
9810a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
98278910aadSHong Zhang   C->assembled           = PETSC_TRUE;
98378910aadSHong Zhang   C->preallocated        = PETSC_TRUE;
98426fbe8dcSKarl Rupp 
9859566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(C->rmap->N));
9863cea4cbeSHong Zhang   if (sctx.nshift) {
987f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
9889566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shiftnz tries %" PetscInt_FMT ", shift_amount %g\n", sctx.nshift, (double)sctx.shift_amount));
989f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
9909566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "number of shiftpd tries %" PetscInt_FMT ", shift_amount %g\n", sctx.nshift, (double)sctx.shift_amount));
99178910aadSHong Zhang     }
99278910aadSHong Zhang   }
993c05c3958SHong Zhang   PetscFunctionReturn(0);
994c05c3958SHong Zhang }
995c05c3958SHong Zhang 
996c6db04a5SJed Brown #include <petscbt.h>
997c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
998d71ae5a4SJacob Faibussowitsch PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact, Mat A, IS perm, const MatFactorInfo *info)
999d71ae5a4SJacob Faibussowitsch {
100078910aadSHong Zhang   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data;
100178910aadSHong Zhang   Mat_SeqSBAIJ      *b;
100278910aadSHong Zhang   Mat                B;
100348dd3d27SHong Zhang   PetscBool          perm_identity, missing;
10045d0c19d7SBarry Smith   PetscInt           reallocs = 0, i, *ai = a->i, *aj = a->j, am = a->mbs, bs = A->rmap->bs, *ui;
10055d0c19d7SBarry Smith   const PetscInt    *rip;
100678910aadSHong Zhang   PetscInt           jmin, jmax, nzk, k, j, *jl, prow, *il, nextprow;
10070298fd71SBarry Smith   PetscInt           nlnk, *lnk, *lnk_lvl = NULL, ncols, ncols_upper, *cols, *cols_lvl, *uj, **uj_ptr, **uj_lvl_ptr;
100878910aadSHong Zhang   PetscReal          fill = info->fill, levels = info->levels;
10090298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL, current_space = NULL;
10100298fd71SBarry Smith   PetscFreeSpaceList free_space_lvl = NULL, current_space_lvl = NULL;
101178910aadSHong Zhang   PetscBT            lnkbt;
101278910aadSHong Zhang 
1013c05c3958SHong Zhang   PetscFunctionBegin;
10149566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal(A, &missing, &i));
101528b400f6SJacob Faibussowitsch   PetscCheck(!missing, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal entry %" PetscInt_FMT, i);
101648dd3d27SHong Zhang 
10176ad2eaddSHong Zhang   if (bs > 1) {
101848a46eb9SPierre Jolivet     if (!a->sbaijMat) PetscCall(MatConvert(A, MATSEQSBAIJ, MAT_INITIAL_MATRIX, &a->sbaijMat));
1019719d5645SBarry Smith     (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */
102026fbe8dcSKarl Rupp 
10219566063dSJacob Faibussowitsch     PetscCall(MatICCFactorSymbolic(fact, a->sbaijMat, perm, info));
10226ad2eaddSHong Zhang     PetscFunctionReturn(0);
10236ad2eaddSHong Zhang   }
10246ad2eaddSHong Zhang 
10259566063dSJacob Faibussowitsch   PetscCall(ISIdentity(perm, &perm_identity));
10269566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(perm, &rip));
102778910aadSHong Zhang 
102878910aadSHong Zhang   /* special case that simply copies fill pattern */
102978910aadSHong Zhang   if (!levels && perm_identity) {
10309566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(am + 1, &ui));
103126fbe8dcSKarl Rupp     for (i = 0; i < am; i++) ui[i] = ai[i + 1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
1032719d5645SBarry Smith     B = fact;
10339566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(B, 1, 0, ui));
103478910aadSHong Zhang 
103578910aadSHong Zhang     b  = (Mat_SeqSBAIJ *)B->data;
103678910aadSHong Zhang     uj = b->j;
103778910aadSHong Zhang     for (i = 0; i < am; i++) {
103878910aadSHong Zhang       aj = a->j + a->diag[i];
103926fbe8dcSKarl Rupp       for (j = 0; j < ui[i]; j++) *uj++ = *aj++;
104078910aadSHong Zhang       b->ilen[i] = ui[i];
104178910aadSHong Zhang     }
10429566063dSJacob Faibussowitsch     PetscCall(PetscFree(ui));
104326fbe8dcSKarl Rupp 
1044d5f3da31SBarry Smith     B->factortype = MAT_FACTOR_NONE;
104526fbe8dcSKarl Rupp 
10469566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
10479566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
1048d5f3da31SBarry Smith     B->factortype = MAT_FACTOR_ICC;
104978910aadSHong Zhang 
105078910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
105178910aadSHong Zhang     PetscFunctionReturn(0);
105278910aadSHong Zhang   }
105378910aadSHong Zhang 
105478910aadSHong Zhang   /* initialization */
10559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(am + 1, &ui));
1056e60cf9a0SBarry Smith   ui[0] = 0;
10579566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * am + 1, &cols_lvl));
105878910aadSHong Zhang 
105978910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
106078910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
10619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(am, &uj_ptr, am, &uj_lvl_ptr, am, &il, am, &jl));
106278910aadSHong Zhang   for (i = 0; i < am; i++) {
10639371c9d4SSatish Balay     jl[i] = am;
10649371c9d4SSatish Balay     il[i] = 0;
106578910aadSHong Zhang   }
106678910aadSHong Zhang 
106778910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
106878910aadSHong Zhang   nlnk = am + 1;
10699566063dSJacob Faibussowitsch   PetscCall(PetscIncompleteLLCreate(am, am, nlnk, lnk, lnk_lvl, lnkbt));
107078910aadSHong Zhang 
107195b5ac22SHong Zhang   /* initial FreeSpace size is fill*(ai[am]+am)/2 */
10729566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill, PetscIntSumTruncate(ai[am] / 2, am / 2)), &free_space));
107326fbe8dcSKarl Rupp 
107478910aadSHong Zhang   current_space = free_space;
107526fbe8dcSKarl Rupp 
10769566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill, PetscIntSumTruncate(ai[am] / 2, am / 2)), &free_space_lvl));
107778910aadSHong Zhang   current_space_lvl = free_space_lvl;
107878910aadSHong Zhang 
107978910aadSHong Zhang   for (k = 0; k < am; k++) { /* for each active row k */
108078910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
108178910aadSHong Zhang     nzk         = 0;
108278910aadSHong Zhang     ncols       = ai[rip[k] + 1] - ai[rip[k]];
108378910aadSHong Zhang     ncols_upper = 0;
108478910aadSHong Zhang     cols        = cols_lvl + am;
108578910aadSHong Zhang     for (j = 0; j < ncols; j++) {
108678910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
108778910aadSHong Zhang       if (i >= k) { /* only take upper triangular entry */
108878910aadSHong Zhang         cols[ncols_upper]     = i;
108978910aadSHong Zhang         cols_lvl[ncols_upper] = -1; /* initialize level for nonzero entries */
109078910aadSHong Zhang         ncols_upper++;
109178910aadSHong Zhang       }
109278910aadSHong Zhang     }
10939566063dSJacob Faibussowitsch     PetscCall(PetscIncompleteLLAdd(ncols_upper, cols, levels, cols_lvl, am, &nlnk, lnk, lnk_lvl, lnkbt));
109478910aadSHong Zhang     nzk += nlnk;
109578910aadSHong Zhang 
109678910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
109778910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
109878910aadSHong Zhang 
109978910aadSHong Zhang     while (prow < k) {
110078910aadSHong Zhang       nextprow = jl[prow];
110178910aadSHong Zhang 
110278910aadSHong Zhang       /* merge prow into k-th row */
110378910aadSHong Zhang       jmin  = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */
110478910aadSHong Zhang       jmax  = ui[prow + 1];
110578910aadSHong Zhang       ncols = jmax - jmin;
110678910aadSHong Zhang       i     = jmin - ui[prow];
110778910aadSHong Zhang       cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
110878910aadSHong Zhang       for (j = 0; j < ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
11099566063dSJacob Faibussowitsch       PetscCall(PetscIncompleteLLAddSorted(ncols, cols, levels, cols_lvl, am, &nlnk, lnk, lnk_lvl, lnkbt));
111078910aadSHong Zhang       nzk += nlnk;
111178910aadSHong Zhang 
111278910aadSHong Zhang       /* update il and jl for prow */
111378910aadSHong Zhang       if (jmin < jmax) {
111478910aadSHong Zhang         il[prow] = jmin;
111526fbe8dcSKarl Rupp 
11169371c9d4SSatish Balay         j        = *cols;
11179371c9d4SSatish Balay         jl[prow] = jl[j];
11189371c9d4SSatish Balay         jl[j]    = prow;
111978910aadSHong Zhang       }
112078910aadSHong Zhang       prow = nextprow;
112178910aadSHong Zhang     }
112278910aadSHong Zhang 
112378910aadSHong Zhang     /* if free space is not available, make more free space */
112478910aadSHong Zhang     if (current_space->local_remaining < nzk) {
112578910aadSHong Zhang       i = am - k + 1;                                                             /* num of unfactored rows */
1126f91af8c7SBarry Smith       i = PetscMin(PetscIntMultTruncate(i, nzk), PetscIntMultTruncate(i, i - 1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
11279566063dSJacob Faibussowitsch       PetscCall(PetscFreeSpaceGet(i, &current_space));
11289566063dSJacob Faibussowitsch       PetscCall(PetscFreeSpaceGet(i, &current_space_lvl));
112978910aadSHong Zhang       reallocs++;
113078910aadSHong Zhang     }
113178910aadSHong Zhang 
113278910aadSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
11339566063dSJacob Faibussowitsch     PetscCall(PetscIncompleteLLClean(am, am, nzk, lnk, lnk_lvl, current_space->array, current_space_lvl->array, lnkbt));
113478910aadSHong Zhang 
113578910aadSHong Zhang     /* add the k-th row into il and jl */
113678910aadSHong Zhang     if (nzk - 1 > 0) {
113778910aadSHong Zhang       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
11389371c9d4SSatish Balay       jl[k] = jl[i];
11399371c9d4SSatish Balay       jl[i] = k;
114078910aadSHong Zhang       il[k] = ui[k] + 1;
114178910aadSHong Zhang     }
114278910aadSHong Zhang     uj_ptr[k]     = current_space->array;
114378910aadSHong Zhang     uj_lvl_ptr[k] = current_space_lvl->array;
114478910aadSHong Zhang 
114578910aadSHong Zhang     current_space->array += nzk;
114678910aadSHong Zhang     current_space->local_used += nzk;
114778910aadSHong Zhang     current_space->local_remaining -= nzk;
114878910aadSHong Zhang 
114978910aadSHong Zhang     current_space_lvl->array += nzk;
115078910aadSHong Zhang     current_space_lvl->local_used += nzk;
115178910aadSHong Zhang     current_space_lvl->local_remaining -= nzk;
115278910aadSHong Zhang 
115378910aadSHong Zhang     ui[k + 1] = ui[k] + nzk;
115478910aadSHong Zhang   }
115578910aadSHong Zhang 
11569566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(perm, &rip));
11579566063dSJacob Faibussowitsch   PetscCall(PetscFree4(uj_ptr, uj_lvl_ptr, il, jl));
11589566063dSJacob Faibussowitsch   PetscCall(PetscFree(cols_lvl));
115978910aadSHong Zhang 
11609263d837SHong Zhang   /* copy free_space into uj and free free_space; set uj in new datastructure; */
11619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ui[am] + 1, &uj));
11629566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceContiguous(&free_space, uj));
11639566063dSJacob Faibussowitsch   PetscCall(PetscIncompleteLLDestroy(lnk, lnkbt));
11649566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceDestroy(free_space_lvl));
116578910aadSHong Zhang 
116678910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1167719d5645SBarry Smith   B = fact;
11689566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(B, 1, MAT_SKIP_ALLOCATION, NULL));
116978910aadSHong Zhang 
117078910aadSHong Zhang   b               = (Mat_SeqSBAIJ *)B->data;
117178910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
1172e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1173e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
117426fbe8dcSKarl Rupp 
11759566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ui[am] + 1, &b->a));
117626fbe8dcSKarl Rupp 
117778910aadSHong Zhang   b->j             = uj;
117878910aadSHong Zhang   b->i             = ui;
1179f4259b30SLisandro Dalcin   b->diag          = NULL;
1180f4259b30SLisandro Dalcin   b->ilen          = NULL;
1181f4259b30SLisandro Dalcin   b->imax          = NULL;
118278910aadSHong Zhang   b->row           = perm;
118378910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
118426fbe8dcSKarl Rupp 
11859566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)perm));
118626fbe8dcSKarl Rupp 
118778910aadSHong Zhang   b->icol = perm;
118826fbe8dcSKarl Rupp 
11899566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)perm));
11909566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(am + 1, &b->solve_work));
119126fbe8dcSKarl Rupp 
119278910aadSHong Zhang   b->maxnz = b->nz = ui[am];
119378910aadSHong Zhang 
119478910aadSHong Zhang   B->info.factor_mallocs   = reallocs;
119578910aadSHong Zhang   B->info.fill_ratio_given = fill;
119675567043SBarry Smith   if (ai[am] != 0.) {
1197*da81f932SPierre Jolivet     /* nonzeros in lower triangular part of A (including diagonals)= (ai[am]+am)/2 */
119895b5ac22SHong Zhang     B->info.fill_ratio_needed = ((PetscReal)2 * ui[am]) / (ai[am] + am);
119978910aadSHong Zhang   } else {
120078910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
120178910aadSHong Zhang   }
12029263d837SHong Zhang #if defined(PETSC_USE_INFO)
12039263d837SHong Zhang   if (ai[am] != 0) {
120495b5ac22SHong Zhang     PetscReal af = B->info.fill_ratio_needed;
12059566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n", reallocs, (double)fill, (double)af));
12069566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Run with -pc_factor_fill %g or use \n", (double)af));
12079566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "PCFactorSetFill(pc,%g) for best performance.\n", (double)af));
12089263d837SHong Zhang   } else {
12099566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Empty matrix\n"));
12109263d837SHong Zhang   }
12119263d837SHong Zhang #endif
121278910aadSHong Zhang   if (perm_identity) {
12130a3c351aSHong Zhang     B->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
12140a3c351aSHong Zhang     B->ops->solvetranspose        = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
121578910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
121678910aadSHong Zhang   } else {
1217719d5645SBarry Smith     (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
121878910aadSHong Zhang   }
1219c05c3958SHong Zhang   PetscFunctionReturn(0);
1220c05c3958SHong Zhang }
1221c05c3958SHong Zhang 
1222d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact, Mat A, IS perm, const MatFactorInfo *info)
1223d71ae5a4SJacob Faibussowitsch {
122478910aadSHong Zhang   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data;
122578910aadSHong Zhang   Mat_SeqSBAIJ      *b;
122678910aadSHong Zhang   Mat                B;
12279186b105SHong Zhang   PetscBool          perm_identity, missing;
122878910aadSHong Zhang   PetscReal          fill = info->fill;
12295d0c19d7SBarry Smith   const PetscInt    *rip;
12305d0c19d7SBarry Smith   PetscInt           i, mbs = a->mbs, bs = A->rmap->bs, *ai = a->i, *aj = a->j, reallocs = 0, prow;
123178910aadSHong Zhang   PetscInt          *jl, jmin, jmax, nzk, *ui, k, j, *il, nextprow;
123278910aadSHong Zhang   PetscInt           nlnk, *lnk, ncols, ncols_upper, *cols, *uj, **ui_ptr, *uj_ptr;
12330298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL, current_space = NULL;
123478910aadSHong Zhang   PetscBT            lnkbt;
123578910aadSHong Zhang 
1236c05c3958SHong Zhang   PetscFunctionBegin;
12376ad2eaddSHong Zhang   if (bs > 1) { /* convert to seqsbaij */
123848a46eb9SPierre Jolivet     if (!a->sbaijMat) PetscCall(MatConvert(A, MATSEQSBAIJ, MAT_INITIAL_MATRIX, &a->sbaijMat));
1239719d5645SBarry Smith     (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */
124026fbe8dcSKarl Rupp 
12419566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorSymbolic(fact, a->sbaijMat, perm, info));
12426ad2eaddSHong Zhang     PetscFunctionReturn(0);
12436ad2eaddSHong Zhang   }
12446ad2eaddSHong Zhang 
12459566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal(A, &missing, &i));
124628b400f6SJacob Faibussowitsch   PetscCheck(!missing, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal entry %" PetscInt_FMT, i);
12479186b105SHong Zhang 
124878910aadSHong Zhang   /* check whether perm is the identity mapping */
12499566063dSJacob Faibussowitsch   PetscCall(ISIdentity(perm, &perm_identity));
125028b400f6SJacob Faibussowitsch   PetscCheck(perm_identity, PETSC_COMM_SELF, PETSC_ERR_SUP, "Matrix reordering is not supported");
12519566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(perm, &rip));
125278910aadSHong Zhang 
125378910aadSHong Zhang   /* initialization */
12549566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs + 1, &ui));
1255e60cf9a0SBarry Smith   ui[0] = 0;
125678910aadSHong Zhang 
125778910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
125878910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
12599566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(mbs, &ui_ptr, mbs, &il, mbs, &jl, mbs, &cols));
126078910aadSHong Zhang   for (i = 0; i < mbs; i++) {
12619371c9d4SSatish Balay     jl[i] = mbs;
12629371c9d4SSatish Balay     il[i] = 0;
126378910aadSHong Zhang   }
126478910aadSHong Zhang 
126578910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
126678910aadSHong Zhang   nlnk = mbs + 1;
12679566063dSJacob Faibussowitsch   PetscCall(PetscLLCreate(mbs, mbs, nlnk, lnk, lnkbt));
126878910aadSHong Zhang 
12696a69fef8SHong Zhang   /* initial FreeSpace size is fill* (ai[mbs]+mbs)/2 */
12709566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill, PetscIntSumTruncate(ai[mbs] / 2, mbs / 2)), &free_space));
127126fbe8dcSKarl Rupp 
127278910aadSHong Zhang   current_space = free_space;
127378910aadSHong Zhang 
127478910aadSHong Zhang   for (k = 0; k < mbs; k++) { /* for each active row k */
127578910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
127678910aadSHong Zhang     nzk         = 0;
127778910aadSHong Zhang     ncols       = ai[rip[k] + 1] - ai[rip[k]];
1278e60cf9a0SBarry Smith     ncols_upper = 0;
127978910aadSHong Zhang     for (j = 0; j < ncols; j++) {
128078910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
128178910aadSHong Zhang       if (i >= k) { /* only take upper triangular entry */
128278910aadSHong Zhang         cols[ncols_upper] = i;
128378910aadSHong Zhang         ncols_upper++;
128478910aadSHong Zhang       }
128578910aadSHong Zhang     }
12869566063dSJacob Faibussowitsch     PetscCall(PetscLLAdd(ncols_upper, cols, mbs, &nlnk, lnk, lnkbt));
128778910aadSHong Zhang     nzk += nlnk;
128878910aadSHong Zhang 
128978910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
129078910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
129178910aadSHong Zhang 
129278910aadSHong Zhang     while (prow < k) {
129378910aadSHong Zhang       nextprow = jl[prow];
129478910aadSHong Zhang       /* merge prow into k-th row */
129578910aadSHong Zhang       jmin   = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
129678910aadSHong Zhang       jmax   = ui[prow + 1];
129778910aadSHong Zhang       ncols  = jmax - jmin;
129878910aadSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
12999566063dSJacob Faibussowitsch       PetscCall(PetscLLAddSorted(ncols, uj_ptr, mbs, &nlnk, lnk, lnkbt));
130078910aadSHong Zhang       nzk += nlnk;
130178910aadSHong Zhang 
130278910aadSHong Zhang       /* update il and jl for prow */
130378910aadSHong Zhang       if (jmin < jmax) {
130478910aadSHong Zhang         il[prow] = jmin;
130526fbe8dcSKarl Rupp         j        = *uj_ptr;
130626fbe8dcSKarl Rupp         jl[prow] = jl[j];
130726fbe8dcSKarl Rupp         jl[j]    = prow;
130878910aadSHong Zhang       }
130978910aadSHong Zhang       prow = nextprow;
131078910aadSHong Zhang     }
131178910aadSHong Zhang 
131278910aadSHong Zhang     /* if free space is not available, make more free space */
131378910aadSHong Zhang     if (current_space->local_remaining < nzk) {
131478910aadSHong Zhang       i = mbs - k + 1;                                                            /* num of unfactored rows */
1315f91af8c7SBarry Smith       i = PetscMin(PetscIntMultTruncate(i, nzk), PetscIntMultTruncate(i, i - 1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
13169566063dSJacob Faibussowitsch       PetscCall(PetscFreeSpaceGet(i, &current_space));
131778910aadSHong Zhang       reallocs++;
131878910aadSHong Zhang     }
131978910aadSHong Zhang 
132078910aadSHong Zhang     /* copy data into free space, then initialize lnk */
13219566063dSJacob Faibussowitsch     PetscCall(PetscLLClean(mbs, mbs, nzk, lnk, current_space->array, lnkbt));
132278910aadSHong Zhang 
132378910aadSHong Zhang     /* add the k-th row into il and jl */
132478910aadSHong Zhang     if (nzk - 1 > 0) {
132578910aadSHong Zhang       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
13269371c9d4SSatish Balay       jl[k] = jl[i];
13279371c9d4SSatish Balay       jl[i] = k;
132878910aadSHong Zhang       il[k] = ui[k] + 1;
132978910aadSHong Zhang     }
133078910aadSHong Zhang     ui_ptr[k] = current_space->array;
133178910aadSHong Zhang     current_space->array += nzk;
133278910aadSHong Zhang     current_space->local_used += nzk;
133378910aadSHong Zhang     current_space->local_remaining -= nzk;
133478910aadSHong Zhang 
133578910aadSHong Zhang     ui[k + 1] = ui[k] + nzk;
133678910aadSHong Zhang   }
133778910aadSHong Zhang 
13389566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(perm, &rip));
13399566063dSJacob Faibussowitsch   PetscCall(PetscFree4(ui_ptr, il, jl, cols));
134078910aadSHong Zhang 
13419263d837SHong Zhang   /* copy free_space into uj and free free_space; set uj in new datastructure; */
13429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ui[mbs] + 1, &uj));
13439566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceContiguous(&free_space, uj));
13449566063dSJacob Faibussowitsch   PetscCall(PetscLLDestroy(lnk, lnkbt));
134578910aadSHong Zhang 
134678910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1347719d5645SBarry Smith   B = fact;
13489566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(B, bs, MAT_SKIP_ALLOCATION, NULL));
134978910aadSHong Zhang 
135078910aadSHong Zhang   b               = (Mat_SeqSBAIJ *)B->data;
135178910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
1352e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1353e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
135426fbe8dcSKarl Rupp 
13559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ui[mbs] + 1, &b->a));
135626fbe8dcSKarl Rupp 
135778910aadSHong Zhang   b->j             = uj;
135878910aadSHong Zhang   b->i             = ui;
1359f4259b30SLisandro Dalcin   b->diag          = NULL;
1360f4259b30SLisandro Dalcin   b->ilen          = NULL;
1361f4259b30SLisandro Dalcin   b->imax          = NULL;
136278910aadSHong Zhang   b->row           = perm;
136378910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
136426fbe8dcSKarl Rupp 
13659566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)perm));
136678910aadSHong Zhang   b->icol = perm;
13679566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)perm));
13689566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs + 1, &b->solve_work));
136978910aadSHong Zhang   b->maxnz = b->nz = ui[mbs];
137078910aadSHong Zhang 
137178910aadSHong Zhang   B->info.factor_mallocs   = reallocs;
137278910aadSHong Zhang   B->info.fill_ratio_given = fill;
137375567043SBarry Smith   if (ai[mbs] != 0.) {
137495b5ac22SHong Zhang     /* nonzeros in lower triangular part of A = (ai[mbs]+mbs)/2 */
137595b5ac22SHong Zhang     B->info.fill_ratio_needed = ((PetscReal)2 * ui[mbs]) / (ai[mbs] + mbs);
137678910aadSHong Zhang   } else {
137778910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
137878910aadSHong Zhang   }
13799263d837SHong Zhang #if defined(PETSC_USE_INFO)
13809263d837SHong Zhang   if (ai[mbs] != 0.) {
13819263d837SHong Zhang     PetscReal af = B->info.fill_ratio_needed;
13829566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n", reallocs, (double)fill, (double)af));
13839566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Run with -pc_factor_fill %g or use \n", (double)af));
13849566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "PCFactorSetFill(pc,%g) for best performance.\n", (double)af));
13859263d837SHong Zhang   } else {
13869566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Empty matrix\n"));
13879263d837SHong Zhang   }
13889263d837SHong Zhang #endif
138978910aadSHong Zhang   if (perm_identity) {
13906ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
139178910aadSHong Zhang   } else {
13926ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
139378910aadSHong Zhang   }
1394c05c3958SHong Zhang   PetscFunctionReturn(0);
1395c05c3958SHong Zhang }
1396c8342467SHong Zhang 
1397d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_SeqBAIJ_N_NaturalOrdering(Mat A, Vec bb, Vec xx)
1398d71ae5a4SJacob Faibussowitsch {
13991a83e813SShri Abhyankar   Mat_SeqBAIJ       *a  = (Mat_SeqBAIJ *)A->data;
14001a83e813SShri Abhyankar   const PetscInt    *ai = a->i, *aj = a->j, *adiag = a->diag, *vi;
14011a83e813SShri Abhyankar   PetscInt           i, k, n                       = a->mbs;
14021a83e813SShri Abhyankar   PetscInt           nz, bs = A->rmap->bs, bs2 = a->bs2;
1403d9ca1df4SBarry Smith   const MatScalar   *aa = a->a, *v;
1404d9ca1df4SBarry Smith   PetscScalar       *x, *s, *t, *ls;
1405d9ca1df4SBarry Smith   const PetscScalar *b;
14061a83e813SShri Abhyankar 
14071a83e813SShri Abhyankar   PetscFunctionBegin;
14089566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
14099566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
14101a83e813SShri Abhyankar   t = a->solve_work;
14111a83e813SShri Abhyankar 
14121a83e813SShri Abhyankar   /* forward solve the lower triangular */
14139566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t, b, bs)); /* copy 1st block of b to t */
14141a83e813SShri Abhyankar 
14151a83e813SShri Abhyankar   for (i = 1; i < n; i++) {
14161a83e813SShri Abhyankar     v  = aa + bs2 * ai[i];
14171a83e813SShri Abhyankar     vi = aj + ai[i];
14181a83e813SShri Abhyankar     nz = ai[i + 1] - ai[i];
14191a83e813SShri Abhyankar     s  = t + bs * i;
14209566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(s, b + bs * i, bs)); /* copy i_th block of b to t */
14211a83e813SShri Abhyankar     for (k = 0; k < nz; k++) {
142296b95a6bSBarry Smith       PetscKernel_v_gets_v_minus_A_times_w(bs, s, v, t + bs * vi[k]);
14231a83e813SShri Abhyankar       v += bs2;
14241a83e813SShri Abhyankar     }
14251a83e813SShri Abhyankar   }
14261a83e813SShri Abhyankar 
14271a83e813SShri Abhyankar   /* backward solve the upper triangular */
14281a83e813SShri Abhyankar   ls = a->solve_work + A->cmap->n;
14291a83e813SShri Abhyankar   for (i = n - 1; i >= 0; i--) {
14301a83e813SShri Abhyankar     v  = aa + bs2 * (adiag[i + 1] + 1);
14311a83e813SShri Abhyankar     vi = aj + adiag[i + 1] + 1;
14321a83e813SShri Abhyankar     nz = adiag[i] - adiag[i + 1] - 1;
14339566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(ls, t + i * bs, bs));
14341a83e813SShri Abhyankar     for (k = 0; k < nz; k++) {
143596b95a6bSBarry Smith       PetscKernel_v_gets_v_minus_A_times_w(bs, ls, v, t + bs * vi[k]);
14361a83e813SShri Abhyankar       v += bs2;
14371a83e813SShri Abhyankar     }
143896b95a6bSBarry Smith     PetscKernel_w_gets_A_times_v(bs, ls, aa + bs2 * adiag[i], t + i * bs); /* *inv(diagonal[i]) */
14399566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(x + i * bs, t + i * bs, bs));
14401a83e813SShri Abhyankar   }
14411a83e813SShri Abhyankar 
14429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
14439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
14449566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * (a->bs2) * (a->nz) - A->rmap->bs * A->cmap->n));
14451a83e813SShri Abhyankar   PetscFunctionReturn(0);
14461a83e813SShri Abhyankar }
14471a83e813SShri Abhyankar 
1448d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_SeqBAIJ_N(Mat A, Vec bb, Vec xx)
1449d71ae5a4SJacob Faibussowitsch {
145035aa4fcfSShri Abhyankar   Mat_SeqBAIJ       *a     = (Mat_SeqBAIJ *)A->data;
145135aa4fcfSShri Abhyankar   IS                 iscol = a->col, isrow = a->row;
145235aa4fcfSShri Abhyankar   const PetscInt    *r, *c, *rout, *cout, *ai = a->i, *aj = a->j, *adiag = a->diag, *vi;
145335aa4fcfSShri Abhyankar   PetscInt           i, m, n = a->mbs;
145435aa4fcfSShri Abhyankar   PetscInt           nz, bs = A->rmap->bs, bs2 = a->bs2;
1455d9ca1df4SBarry Smith   const MatScalar   *aa = a->a, *v;
1456d9ca1df4SBarry Smith   PetscScalar       *x, *s, *t, *ls;
1457d9ca1df4SBarry Smith   const PetscScalar *b;
145835aa4fcfSShri Abhyankar 
145935aa4fcfSShri Abhyankar   PetscFunctionBegin;
14609566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
14619566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
146235aa4fcfSShri Abhyankar   t = a->solve_work;
146335aa4fcfSShri Abhyankar 
14649371c9d4SSatish Balay   PetscCall(ISGetIndices(isrow, &rout));
14659371c9d4SSatish Balay   r = rout;
14669371c9d4SSatish Balay   PetscCall(ISGetIndices(iscol, &cout));
14679371c9d4SSatish Balay   c = cout;
146835aa4fcfSShri Abhyankar 
146935aa4fcfSShri Abhyankar   /* forward solve the lower triangular */
14709566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t, b + bs * r[0], bs));
147135aa4fcfSShri Abhyankar   for (i = 1; i < n; i++) {
147235aa4fcfSShri Abhyankar     v  = aa + bs2 * ai[i];
147335aa4fcfSShri Abhyankar     vi = aj + ai[i];
147435aa4fcfSShri Abhyankar     nz = ai[i + 1] - ai[i];
147535aa4fcfSShri Abhyankar     s  = t + bs * i;
14769566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(s, b + bs * r[i], bs));
147735aa4fcfSShri Abhyankar     for (m = 0; m < nz; m++) {
147896b95a6bSBarry Smith       PetscKernel_v_gets_v_minus_A_times_w(bs, s, v, t + bs * vi[m]);
147935aa4fcfSShri Abhyankar       v += bs2;
148035aa4fcfSShri Abhyankar     }
148135aa4fcfSShri Abhyankar   }
148235aa4fcfSShri Abhyankar 
148335aa4fcfSShri Abhyankar   /* backward solve the upper triangular */
148435aa4fcfSShri Abhyankar   ls = a->solve_work + A->cmap->n;
148535aa4fcfSShri Abhyankar   for (i = n - 1; i >= 0; i--) {
148635aa4fcfSShri Abhyankar     v  = aa + bs2 * (adiag[i + 1] + 1);
148735aa4fcfSShri Abhyankar     vi = aj + adiag[i + 1] + 1;
148835aa4fcfSShri Abhyankar     nz = adiag[i] - adiag[i + 1] - 1;
14899566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(ls, t + i * bs, bs));
149035aa4fcfSShri Abhyankar     for (m = 0; m < nz; m++) {
149196b95a6bSBarry Smith       PetscKernel_v_gets_v_minus_A_times_w(bs, ls, v, t + bs * vi[m]);
149235aa4fcfSShri Abhyankar       v += bs2;
149335aa4fcfSShri Abhyankar     }
149496b95a6bSBarry Smith     PetscKernel_w_gets_A_times_v(bs, ls, v, t + i * bs); /* *inv(diagonal[i]) */
14959566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(x + bs * c[i], t + i * bs, bs));
149635aa4fcfSShri Abhyankar   }
14979566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &rout));
14989566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscol, &cout));
14999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
15009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
15019566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * (a->bs2) * (a->nz) - A->rmap->bs * A->cmap->n));
150235aa4fcfSShri Abhyankar   PetscFunctionReturn(0);
150335aa4fcfSShri Abhyankar }
150435aa4fcfSShri Abhyankar 
1505a25532f0SBarry Smith /*
1506a25532f0SBarry Smith     For each block in an block array saves the largest absolute value in the block into another array
1507a25532f0SBarry Smith */
1508d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatBlockAbs_private(PetscInt nbs, PetscInt bs2, PetscScalar *blockarray, PetscReal *absarray)
1509d71ae5a4SJacob Faibussowitsch {
15102efa7f71SHong Zhang   PetscInt i, j;
15115fd66863SKarl Rupp 
15122efa7f71SHong Zhang   PetscFunctionBegin;
15139566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(absarray, nbs + 1));
15142efa7f71SHong Zhang   for (i = 0; i < nbs; i++) {
15152efa7f71SHong Zhang     for (j = 0; j < bs2; j++) {
15162efa7f71SHong Zhang       if (absarray[i] < PetscAbsScalar(blockarray[i * nbs + j])) absarray[i] = PetscAbsScalar(blockarray[i * nbs + j]);
15172efa7f71SHong Zhang     }
15182efa7f71SHong Zhang   }
15192efa7f71SHong Zhang   PetscFunctionReturn(0);
15202efa7f71SHong Zhang }
15212efa7f71SHong Zhang 
1522fe97e370SBarry Smith /*
1523fe97e370SBarry Smith      This needs to be renamed and called by the regular MatILUFactor_SeqBAIJ when drop tolerance is used
1524fe97e370SBarry Smith */
1525d71ae5a4SJacob Faibussowitsch PetscErrorCode MatILUDTFactor_SeqBAIJ(Mat A, IS isrow, IS iscol, const MatFactorInfo *info, Mat *fact)
1526d71ae5a4SJacob Faibussowitsch {
15272efa7f71SHong Zhang   Mat             B = *fact;
15282efa7f71SHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b;
15292efa7f71SHong Zhang   IS              isicol;
15302efa7f71SHong Zhang   const PetscInt *r, *ic;
15312efa7f71SHong Zhang   PetscInt        i, mbs = a->mbs, bs = A->rmap->bs, bs2 = a->bs2, *ai = a->i, *aj = a->j, *ajtmp, *adiag;
15322efa7f71SHong Zhang   PetscInt       *bi, *bj, *bdiag;
15332efa7f71SHong Zhang 
15342efa7f71SHong Zhang   PetscInt   row, nzi, nzi_bl, nzi_bu, *im, dtcount, nzi_al, nzi_au;
15352efa7f71SHong Zhang   PetscInt   nlnk, *lnk;
15362efa7f71SHong Zhang   PetscBT    lnkbt;
1537ace3abfcSBarry Smith   PetscBool  row_identity, icol_identity;
15382efa7f71SHong Zhang   MatScalar *aatmp, *pv, *batmp, *ba, *rtmp, *pc, *multiplier, *vtmp;
15392efa7f71SHong Zhang   PetscInt   j, nz, *pj, *bjtmp, k, ncut, *jtmp;
15402efa7f71SHong Zhang 
1541182b8fbaSHong Zhang   PetscReal  dt = info->dt; /* shift=info->shiftamount; */
15422efa7f71SHong Zhang   PetscInt   nnz_max;
1543ace3abfcSBarry Smith   PetscBool  missing;
1544021822bcSHong Zhang   PetscReal *vtmp_abs;
154597ef94ebSSatish Balay   MatScalar *v_work;
154697ef94ebSSatish Balay   PetscInt  *v_pivots;
15475f8bbccaSHong Zhang   PetscBool  allowzeropivot, zeropivotdetected = PETSC_FALSE;
15482efa7f71SHong Zhang 
1549c8342467SHong Zhang   PetscFunctionBegin;
15502efa7f71SHong Zhang   /* ------- symbolic factorization, can be reused ---------*/
15519566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal(A, &missing, &i));
155228b400f6SJacob Faibussowitsch   PetscCheck(!missing, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal entry %" PetscInt_FMT, i);
15532efa7f71SHong Zhang   adiag = a->diag;
15542efa7f71SHong Zhang 
15559566063dSJacob Faibussowitsch   PetscCall(ISInvertPermutation(iscol, PETSC_DECIDE, &isicol));
15562efa7f71SHong Zhang 
15572efa7f71SHong Zhang   /* bdiag is location of diagonal in factor */
15589566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs + 1, &bdiag));
15592efa7f71SHong Zhang 
15602efa7f71SHong Zhang   /* allocate row pointers bi */
15619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * mbs + 2, &bi));
15622efa7f71SHong Zhang 
15632efa7f71SHong Zhang   /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */
15642efa7f71SHong Zhang   dtcount = (PetscInt)info->dtcount;
15656bce7ff8SHong Zhang   if (dtcount > mbs - 1) dtcount = mbs - 1;
15666bce7ff8SHong Zhang   nnz_max = ai[mbs] + 2 * mbs * dtcount + 2;
15676da515a1SHong Zhang   /* printf("MatILUDTFactor_SeqBAIJ, bs %d, ai[mbs] %d, nnz_max  %d, dtcount %d\n",bs,ai[mbs],nnz_max,dtcount); */
15689566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nnz_max, &bj));
15696bce7ff8SHong Zhang   nnz_max = nnz_max * bs2;
15709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nnz_max, &ba));
15712efa7f71SHong Zhang 
15722efa7f71SHong Zhang   /* put together the new matrix */
15739566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(B, bs, MAT_SKIP_ALLOCATION, NULL));
157426fbe8dcSKarl Rupp 
15752efa7f71SHong Zhang   b               = (Mat_SeqBAIJ *)(B)->data;
15762efa7f71SHong Zhang   b->free_a       = PETSC_TRUE;
15772efa7f71SHong Zhang   b->free_ij      = PETSC_TRUE;
15782efa7f71SHong Zhang   b->singlemalloc = PETSC_FALSE;
157926fbe8dcSKarl Rupp 
15802efa7f71SHong Zhang   b->a    = ba;
15812efa7f71SHong Zhang   b->j    = bj;
15822efa7f71SHong Zhang   b->i    = bi;
15832efa7f71SHong Zhang   b->diag = bdiag;
1584f4259b30SLisandro Dalcin   b->ilen = NULL;
1585f4259b30SLisandro Dalcin   b->imax = NULL;
15862efa7f71SHong Zhang   b->row  = isrow;
15872efa7f71SHong Zhang   b->col  = iscol;
158826fbe8dcSKarl Rupp 
15899566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)isrow));
15909566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)iscol));
159126fbe8dcSKarl Rupp 
15922efa7f71SHong Zhang   b->icol = isicol;
15939566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs * (mbs + 1), &b->solve_work));
15942efa7f71SHong Zhang   b->maxnz = nnz_max / bs2;
15952efa7f71SHong Zhang 
1596d5f3da31SBarry Smith   (B)->factortype            = MAT_FACTOR_ILUDT;
15972efa7f71SHong Zhang   (B)->info.factor_mallocs   = 0;
15982efa7f71SHong Zhang   (B)->info.fill_ratio_given = ((PetscReal)nnz_max) / ((PetscReal)(ai[mbs] * bs2));
15992efa7f71SHong Zhang   /* ------- end of symbolic factorization ---------*/
16009566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
16019566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
16022efa7f71SHong Zhang 
16032efa7f71SHong Zhang   /* linked list for storing column indices of the active row */
16042efa7f71SHong Zhang   nlnk = mbs + 1;
16059566063dSJacob Faibussowitsch   PetscCall(PetscLLCreate(mbs, mbs, nlnk, lnk, lnkbt));
16062efa7f71SHong Zhang 
16072efa7f71SHong Zhang   /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */
16089566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(mbs, &im, mbs, &jtmp));
16092efa7f71SHong Zhang   /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */
16109566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(mbs * bs2, &rtmp, mbs * bs2, &vtmp));
16119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs + 1, &vtmp_abs));
16129566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(bs, &v_work, bs2, &multiplier, bs, &v_pivots));
16132efa7f71SHong Zhang 
16145f8bbccaSHong Zhang   allowzeropivot  = PetscNot(A->erroriffailure);
1615e60cf9a0SBarry Smith   bi[0]           = 0;
16162efa7f71SHong Zhang   bdiag[0]        = (nnz_max / bs2) - 1; /* location of diagonal in factor B */
16176bce7ff8SHong Zhang   bi[2 * mbs + 1] = bdiag[0] + 1;        /* endof bj and ba array */
16182efa7f71SHong Zhang   for (i = 0; i < mbs; i++) {
16192efa7f71SHong Zhang     /* copy initial fill into linked list */
16202efa7f71SHong Zhang     nzi = ai[r[i] + 1] - ai[r[i]];
162128b400f6SJacob 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);
16222efa7f71SHong Zhang     nzi_al = adiag[r[i]] - ai[r[i]];
16232efa7f71SHong Zhang     nzi_au = ai[r[i] + 1] - adiag[r[i]] - 1;
16242efa7f71SHong Zhang 
16252efa7f71SHong Zhang     /* load in initial unfactored row */
16262efa7f71SHong Zhang     ajtmp = aj + ai[r[i]];
16279566063dSJacob Faibussowitsch     PetscCall(PetscLLAddPerm(nzi, ajtmp, ic, mbs, &nlnk, lnk, lnkbt));
16289566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(rtmp, mbs * bs2));
16292efa7f71SHong Zhang     aatmp = a->a + bs2 * ai[r[i]];
16309566063dSJacob Faibussowitsch     for (j = 0; j < nzi; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ic[ajtmp[j]], aatmp + bs2 * j, bs2));
16312efa7f71SHong Zhang 
16322efa7f71SHong Zhang     /* add pivot rows into linked list */
16332efa7f71SHong Zhang     row = lnk[mbs];
16342efa7f71SHong Zhang     while (row < i) {
16352efa7f71SHong Zhang       nzi_bl = bi[row + 1] - bi[row] + 1;
16362efa7f71SHong Zhang       bjtmp  = bj + bdiag[row + 1] + 1; /* points to 1st column next to the diagonal in U */
16379566063dSJacob Faibussowitsch       PetscCall(PetscLLAddSortedLU(bjtmp, row, &nlnk, lnk, lnkbt, i, nzi_bl, im));
16382efa7f71SHong Zhang       nzi += nlnk;
16392efa7f71SHong Zhang       row = lnk[row];
16402efa7f71SHong Zhang     }
16412efa7f71SHong Zhang 
16422efa7f71SHong Zhang     /* copy data from lnk into jtmp, then initialize lnk */
16439566063dSJacob Faibussowitsch     PetscCall(PetscLLClean(mbs, mbs, nzi, lnk, jtmp, lnkbt));
16442efa7f71SHong Zhang 
16452efa7f71SHong Zhang     /* numerical factorization */
16462efa7f71SHong Zhang     bjtmp = jtmp;
16472efa7f71SHong Zhang     row   = *bjtmp++; /* 1st pivot row */
16482efa7f71SHong Zhang 
16492efa7f71SHong Zhang     while (row < i) {
16502efa7f71SHong Zhang       pc = rtmp + bs2 * row;
16512efa7f71SHong Zhang       pv = ba + bs2 * bdiag[row];                           /* inv(diag) of the pivot row */
165296b95a6bSBarry Smith       PetscKernel_A_gets_A_times_B(bs, pc, pv, multiplier); /* pc= multiplier = pc*inv(diag[row]) */
16539566063dSJacob Faibussowitsch       PetscCall(MatBlockAbs_private(1, bs2, pc, vtmp_abs));
16542efa7f71SHong Zhang       if (vtmp_abs[0] > dt) {         /* apply tolerance dropping rule */
16552efa7f71SHong Zhang         pj = bj + bdiag[row + 1] + 1; /* point to 1st entry of U(row,:) */
16562efa7f71SHong Zhang         pv = ba + bs2 * (bdiag[row + 1] + 1);
16572efa7f71SHong Zhang         nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries in U(row,:), excluding diagonal */
1658ad540459SPierre Jolivet         for (j = 0; j < nz; j++) PetscKernel_A_gets_A_minus_B_times_C(bs, rtmp + bs2 * pj[j], pc, pv + bs2 * j);
16599566063dSJacob Faibussowitsch         /* PetscCall(PetscLogFlops(bslog*(nz+1.0)-bs)); */
16602efa7f71SHong Zhang       }
16612efa7f71SHong Zhang       row = *bjtmp++;
16622efa7f71SHong Zhang     }
16632efa7f71SHong Zhang 
16642efa7f71SHong Zhang     /* copy sparse rtmp into contiguous vtmp; separate L and U part */
16659371c9d4SSatish Balay     nzi_bl = 0;
16669371c9d4SSatish Balay     j      = 0;
16672efa7f71SHong Zhang     while (jtmp[j] < i) { /* L-part. Note: jtmp is sorted */
16689566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(vtmp + bs2 * j, rtmp + bs2 * jtmp[j], bs2));
16699371c9d4SSatish Balay       nzi_bl++;
16709371c9d4SSatish Balay       j++;
16712efa7f71SHong Zhang     }
16722efa7f71SHong Zhang     nzi_bu = nzi - nzi_bl - 1;
16732efa7f71SHong Zhang 
16742efa7f71SHong Zhang     while (j < nzi) { /* U-part */
16759566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(vtmp + bs2 * j, rtmp + bs2 * jtmp[j], bs2));
16762efa7f71SHong Zhang       j++;
16772efa7f71SHong Zhang     }
16782efa7f71SHong Zhang 
16799566063dSJacob Faibussowitsch     PetscCall(MatBlockAbs_private(nzi, bs2, vtmp, vtmp_abs));
16805f8bbccaSHong Zhang 
16812efa7f71SHong Zhang     bjtmp = bj + bi[i];
16822efa7f71SHong Zhang     batmp = ba + bs2 * bi[i];
16832efa7f71SHong Zhang     /* apply level dropping rule to L part */
16842efa7f71SHong Zhang     ncut = nzi_al + dtcount;
16852efa7f71SHong Zhang     if (ncut < nzi_bl) {
16869566063dSJacob Faibussowitsch       PetscCall(PetscSortSplitReal(ncut, nzi_bl, vtmp_abs, jtmp));
16879566063dSJacob Faibussowitsch       PetscCall(PetscSortIntWithScalarArray(ncut, jtmp, vtmp));
16882efa7f71SHong Zhang     } else {
16892efa7f71SHong Zhang       ncut = nzi_bl;
16902efa7f71SHong Zhang     }
16912efa7f71SHong Zhang     for (j = 0; j < ncut; j++) {
16922efa7f71SHong Zhang       bjtmp[j] = jtmp[j];
16939566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(batmp + bs2 * j, rtmp + bs2 * bjtmp[j], bs2));
16942efa7f71SHong Zhang     }
16952efa7f71SHong Zhang     bi[i + 1] = bi[i] + ncut;
16962efa7f71SHong Zhang     nzi       = ncut + 1;
16972efa7f71SHong Zhang 
16982efa7f71SHong Zhang     /* apply level dropping rule to U part */
16992efa7f71SHong Zhang     ncut = nzi_au + dtcount;
17002efa7f71SHong Zhang     if (ncut < nzi_bu) {
17019566063dSJacob Faibussowitsch       PetscCall(PetscSortSplitReal(ncut, nzi_bu, vtmp_abs + nzi_bl + 1, jtmp + nzi_bl + 1));
17029566063dSJacob Faibussowitsch       PetscCall(PetscSortIntWithScalarArray(ncut, jtmp + nzi_bl + 1, vtmp + nzi_bl + 1));
17032efa7f71SHong Zhang     } else {
17042efa7f71SHong Zhang       ncut = nzi_bu;
17052efa7f71SHong Zhang     }
17062efa7f71SHong Zhang     nzi += ncut;
17072efa7f71SHong Zhang 
17082efa7f71SHong Zhang     /* mark bdiagonal */
17092efa7f71SHong Zhang     bdiag[i + 1]    = bdiag[i] - (ncut + 1);
17106bce7ff8SHong Zhang     bi[2 * mbs - i] = bi[2 * mbs - i + 1] - (ncut + 1);
17116bce7ff8SHong Zhang 
17122efa7f71SHong Zhang     bjtmp = bj + bdiag[i];
17132efa7f71SHong Zhang     batmp = ba + bs2 * bdiag[i];
17149566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(batmp, rtmp + bs2 * i, bs2));
17152efa7f71SHong Zhang     *bjtmp = i;
17165f8bbccaSHong Zhang 
17172efa7f71SHong Zhang     bjtmp = bj + bdiag[i + 1] + 1;
17182efa7f71SHong Zhang     batmp = ba + (bdiag[i + 1] + 1) * bs2;
17192efa7f71SHong Zhang 
17202efa7f71SHong Zhang     for (k = 0; k < ncut; k++) {
17212efa7f71SHong Zhang       bjtmp[k] = jtmp[nzi_bl + 1 + k];
17229566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(batmp + bs2 * k, rtmp + bs2 * bjtmp[k], bs2));
17232efa7f71SHong Zhang     }
17242efa7f71SHong Zhang 
17252efa7f71SHong Zhang     im[i] = nzi; /* used by PetscLLAddSortedLU() */
17262efa7f71SHong Zhang 
1727a5b23f4aSJose E. Roman     /* invert diagonal block for simpler triangular solves - add shift??? */
17282efa7f71SHong Zhang     batmp = ba + bs2 * bdiag[i];
17295f8bbccaSHong Zhang 
17309566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A(bs, batmp, v_pivots, v_work, allowzeropivot, &zeropivotdetected));
17317b6c816cSBarry Smith     if (zeropivotdetected) B->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
17322efa7f71SHong Zhang   } /* for (i=0; i<mbs; i++) */
17339566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v_work, multiplier, v_pivots));
17342efa7f71SHong Zhang 
17356da515a1SHong Zhang   /* printf("end of L %d, beginning of U %d\n",bi[mbs],bdiag[mbs]); */
173694bad497SJacob 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]);
17372efa7f71SHong Zhang 
17389566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
17399566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
17402efa7f71SHong Zhang 
17419566063dSJacob Faibussowitsch   PetscCall(PetscLLDestroy(lnk, lnkbt));
17422efa7f71SHong Zhang 
17439566063dSJacob Faibussowitsch   PetscCall(PetscFree2(im, jtmp));
17449566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rtmp, vtmp));
17452efa7f71SHong Zhang 
17469566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(bs2 * B->cmap->n));
17472efa7f71SHong Zhang   b->maxnz = b->nz = bi[mbs] + bdiag[0] - bdiag[mbs];
17482efa7f71SHong Zhang 
17499566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isrow, &row_identity));
17509566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isicol, &icol_identity));
17512efa7f71SHong Zhang   if (row_identity && icol_identity) {
17524dd39f65SShri Abhyankar     B->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering;
17532efa7f71SHong Zhang   } else {
17544dd39f65SShri Abhyankar     B->ops->solve = MatSolve_SeqBAIJ_N;
17552efa7f71SHong Zhang   }
17562efa7f71SHong Zhang 
1757f4259b30SLisandro Dalcin   B->ops->solveadd          = NULL;
1758f4259b30SLisandro Dalcin   B->ops->solvetranspose    = NULL;
1759f4259b30SLisandro Dalcin   B->ops->solvetransposeadd = NULL;
1760f4259b30SLisandro Dalcin   B->ops->matsolve          = NULL;
17612efa7f71SHong Zhang   B->assembled              = PETSC_TRUE;
17622efa7f71SHong Zhang   B->preallocated           = PETSC_TRUE;
1763c8342467SHong Zhang   PetscFunctionReturn(0);
1764c8342467SHong Zhang }
1765