xref: /petsc/src/mat/impls/baij/seq/baijfact.c (revision b94d7ded0a05f1bbd5e48daa6f92b28259c75b44)
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 
84dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat B,Mat A,const MatFactorInfo *info)
9b588c5a2SHong Zhang {
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];
37b588c5a2SHong Zhang     for  (j=0; j<nz; j++) {
389566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(rtmp+bs2*bjtmp[j],bs2));
39b588c5a2SHong Zhang     }
40b588c5a2SHong Zhang 
41b588c5a2SHong Zhang     /* U part */
420c4413a7SShri Abhyankar     nz    = bdiag[i] - bdiag[i+1];
430c4413a7SShri Abhyankar     bjtmp = bj + bdiag[i+1]+1;
440c4413a7SShri Abhyankar     for  (j=0; j<nz; j++) {
459566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(rtmp+bs2*bjtmp[j],bs2));
460c4413a7SShri Abhyankar     }
470c4413a7SShri Abhyankar 
480c4413a7SShri Abhyankar     /* load in initial (unfactored row) */
490c4413a7SShri Abhyankar     nz    = ai[r[i]+1] - ai[r[i]];
500c4413a7SShri Abhyankar     ajtmp = aj + ai[r[i]];
510c4413a7SShri Abhyankar     v     = aa + bs2*ai[r[i]];
520c4413a7SShri Abhyankar     for (j=0; j<nz; j++) {
539566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2));
540c4413a7SShri Abhyankar     }
550c4413a7SShri Abhyankar 
560c4413a7SShri Abhyankar     /* elimination */
570c4413a7SShri Abhyankar     bjtmp = bj + bi[i];
580c4413a7SShri Abhyankar     nzL   = bi[i+1] - bi[i];
590c4413a7SShri Abhyankar     for (k=0; k < nzL; k++) {
600c4413a7SShri Abhyankar       row = bjtmp[k];
610c4413a7SShri Abhyankar       pc  = rtmp + bs2*row;
62c35f09e5SBarry Smith       for (flg=0,j=0; j<bs2; j++) {
63c35f09e5SBarry Smith         if (pc[j] != (PetscScalar)0.0) {
64c35f09e5SBarry Smith           flg = 1;
65c35f09e5SBarry Smith           break;
66c35f09e5SBarry Smith         }
67c35f09e5SBarry Smith       }
680c4413a7SShri Abhyankar       if (flg) {
690c4413a7SShri Abhyankar         pv = b->a + bs2*bdiag[row];
7096b95a6bSBarry Smith         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
719566063dSJacob Faibussowitsch         PetscCall(PetscKernel_A_gets_A_times_B_2(pc,pv,mwork));
720c4413a7SShri Abhyankar 
73a5b23f4aSJose E. Roman         pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
740c4413a7SShri Abhyankar         pv = b->a + bs2*(bdiag[row+1]+1);
750c4413a7SShri Abhyankar         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
760c4413a7SShri Abhyankar         for (j=0; j<nz; j++) {
7796b95a6bSBarry Smith           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
780c4413a7SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
790c4413a7SShri Abhyankar           v    = rtmp + 4*pj[j];
809566063dSJacob Faibussowitsch           PetscCall(PetscKernel_A_gets_A_minus_B_times_C_2(v,pc,pv));
810c4413a7SShri Abhyankar           pv  += 4;
820c4413a7SShri Abhyankar         }
839566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(16.0*nz+12)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
840c4413a7SShri Abhyankar       }
850c4413a7SShri Abhyankar     }
860c4413a7SShri Abhyankar 
870c4413a7SShri Abhyankar     /* finished row so stick it into b->a */
880c4413a7SShri Abhyankar     /* L part */
890c4413a7SShri Abhyankar     pv = b->a + bs2*bi[i];
900c4413a7SShri Abhyankar     pj = b->j + bi[i];
910c4413a7SShri Abhyankar     nz = bi[i+1] - bi[i];
920c4413a7SShri Abhyankar     for (j=0; j<nz; j++) {
939566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2));
940c4413a7SShri Abhyankar     }
950c4413a7SShri Abhyankar 
96a5b23f4aSJose E. Roman     /* Mark diagonal and invert diagonal for simpler triangular solves */
970c4413a7SShri Abhyankar     pv   = b->a + bs2*bdiag[i];
980c4413a7SShri Abhyankar     pj   = b->j + bdiag[i];
999566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(pv,rtmp+bs2*pj[0],bs2));
1009566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_2(pv,shift,allowzeropivot,&zeropivotdetected));
1017b6c816cSBarry Smith     if (zeropivotdetected) B->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1020c4413a7SShri Abhyankar 
1030c4413a7SShri Abhyankar     /* U part */
1040c4413a7SShri Abhyankar     pv = b->a + bs2*(bdiag[i+1]+1);
1050c4413a7SShri Abhyankar     pj = b->j + bdiag[i+1]+1;
1060c4413a7SShri Abhyankar     nz = bdiag[i] - bdiag[i+1] - 1;
1070c4413a7SShri Abhyankar     for (j=0; j<nz; j++) {
1089566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2));
1090c4413a7SShri Abhyankar     }
1100c4413a7SShri Abhyankar   }
1110c4413a7SShri Abhyankar 
1129566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rtmp,mwork));
1139566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol,&ic));
1149566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow,&r));
11526fbe8dcSKarl Rupp 
1164dd39f65SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_2;
1174dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2;
1180c4413a7SShri Abhyankar   C->assembled           = PETSC_TRUE;
11926fbe8dcSKarl Rupp 
1209566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333*2*2*2*n)); /* from inverting diagonal blocks */
1210c4413a7SShri Abhyankar   PetscFunctionReturn(0);
1220c4413a7SShri Abhyankar }
1230c4413a7SShri Abhyankar 
1244dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
1254c000e2eSHong Zhang {
1264c000e2eSHong Zhang   Mat            C =B;
1274c000e2eSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
128bbd65245SShri Abhyankar   PetscInt       i,j,k,nz,nzL,row,*pj;
129bbd65245SShri Abhyankar   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,bs2=a->bs2;
130bbd65245SShri Abhyankar   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag;
131bbd65245SShri Abhyankar   MatScalar      *rtmp,*pc,*mwork,*pv;
132bbd65245SShri Abhyankar   MatScalar      *aa=a->a,*v;
133bbd65245SShri Abhyankar   PetscInt       flg;
134182b8fbaSHong Zhang   PetscReal      shift = info->shiftamount;
135a455e926SHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
1364c000e2eSHong Zhang 
1374c000e2eSHong Zhang   PetscFunctionBegin;
1380164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
1390164db54SHong Zhang 
1404c000e2eSHong Zhang   /* generate work space needed by the factorization */
1419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(bs2*n,&rtmp,bs2,&mwork));
1429566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(rtmp,bs2*n));
1434c000e2eSHong Zhang 
1444c000e2eSHong Zhang   for (i=0; i<n; i++) {
1454c000e2eSHong Zhang     /* zero rtmp */
1464c000e2eSHong Zhang     /* L part */
1474c000e2eSHong Zhang     nz    = bi[i+1] - bi[i];
1484c000e2eSHong Zhang     bjtmp = bj + bi[i];
1494c000e2eSHong Zhang     for  (j=0; j<nz; j++) {
1509566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(rtmp+bs2*bjtmp[j],bs2));
1514c000e2eSHong Zhang     }
1524c000e2eSHong Zhang 
1534c000e2eSHong Zhang     /* U part */
154b2b2dd24SShri Abhyankar     nz    = bdiag[i] - bdiag[i+1];
155b2b2dd24SShri Abhyankar     bjtmp = bj + bdiag[i+1]+1;
156b2b2dd24SShri Abhyankar     for  (j=0; j<nz; j++) {
1579566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(rtmp+bs2*bjtmp[j],bs2));
158b2b2dd24SShri Abhyankar     }
159b2b2dd24SShri Abhyankar 
160b2b2dd24SShri Abhyankar     /* load in initial (unfactored row) */
161b2b2dd24SShri Abhyankar     nz    = ai[i+1] - ai[i];
162b2b2dd24SShri Abhyankar     ajtmp = aj + ai[i];
163b2b2dd24SShri Abhyankar     v     = aa + bs2*ai[i];
164b2b2dd24SShri Abhyankar     for (j=0; j<nz; j++) {
1659566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2));
166b2b2dd24SShri Abhyankar     }
167b2b2dd24SShri Abhyankar 
168b2b2dd24SShri Abhyankar     /* elimination */
169b2b2dd24SShri Abhyankar     bjtmp = bj + bi[i];
170b2b2dd24SShri Abhyankar     nzL   = bi[i+1] - bi[i];
171b2b2dd24SShri Abhyankar     for (k=0; k < nzL; k++) {
172b2b2dd24SShri Abhyankar       row = bjtmp[k];
173b2b2dd24SShri Abhyankar       pc  = rtmp + bs2*row;
174c35f09e5SBarry Smith       for (flg=0,j=0; j<bs2; j++) {
175c35f09e5SBarry Smith         if (pc[j]!=(PetscScalar)0.0) {
176c35f09e5SBarry Smith           flg = 1;
177c35f09e5SBarry Smith           break;
178c35f09e5SBarry Smith         }
179c35f09e5SBarry Smith       }
180b2b2dd24SShri Abhyankar       if (flg) {
181b2b2dd24SShri Abhyankar         pv = b->a + bs2*bdiag[row];
18296b95a6bSBarry Smith         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
1839566063dSJacob Faibussowitsch         PetscCall(PetscKernel_A_gets_A_times_B_2(pc,pv,mwork));
184b2b2dd24SShri Abhyankar 
185b2b2dd24SShri Abhyankar         pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
186b2b2dd24SShri Abhyankar         pv = b->a + bs2*(bdiag[row+1]+1);
187b2b2dd24SShri Abhyankar         nz = bdiag[row]-bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */
188b2b2dd24SShri Abhyankar         for (j=0; j<nz; j++) {
18996b95a6bSBarry Smith           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
190b2b2dd24SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
191b2b2dd24SShri Abhyankar           v    = rtmp + 4*pj[j];
1929566063dSJacob Faibussowitsch           PetscCall(PetscKernel_A_gets_A_minus_B_times_C_2(v,pc,pv));
193b2b2dd24SShri Abhyankar           pv  += 4;
194b2b2dd24SShri Abhyankar         }
1959566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(16.0*nz+12)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
196b2b2dd24SShri Abhyankar       }
197b2b2dd24SShri Abhyankar     }
198b2b2dd24SShri Abhyankar 
199b2b2dd24SShri Abhyankar     /* finished row so stick it into b->a */
200b2b2dd24SShri Abhyankar     /* L part */
201b2b2dd24SShri Abhyankar     pv = b->a + bs2*bi[i];
202b2b2dd24SShri Abhyankar     pj = b->j + bi[i];
203b2b2dd24SShri Abhyankar     nz = bi[i+1] - bi[i];
204b2b2dd24SShri Abhyankar     for (j=0; j<nz; j++) {
2059566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2));
206b2b2dd24SShri Abhyankar     }
207b2b2dd24SShri Abhyankar 
208a5b23f4aSJose E. Roman     /* Mark diagonal and invert diagonal for simpler triangular solves */
209b2b2dd24SShri Abhyankar     pv   = b->a + bs2*bdiag[i];
210b2b2dd24SShri Abhyankar     pj   = b->j + bdiag[i];
2119566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(pv,rtmp+bs2*pj[0],bs2));
2129566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_2(pv,shift,allowzeropivot,&zeropivotdetected));
2137b6c816cSBarry Smith     if (zeropivotdetected) B->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
214b2b2dd24SShri Abhyankar 
215b2b2dd24SShri Abhyankar     /* U part */
216b2b2dd24SShri Abhyankar     /*
217b2b2dd24SShri Abhyankar     pv = b->a + bs2*bi[2*n-i];
218b2b2dd24SShri Abhyankar     pj = b->j + bi[2*n-i];
219b2b2dd24SShri Abhyankar     nz = bi[2*n-i+1] - bi[2*n-i] - 1;
220b2b2dd24SShri Abhyankar     */
221b2b2dd24SShri Abhyankar     pv = b->a + bs2*(bdiag[i+1]+1);
222b2b2dd24SShri Abhyankar     pj = b->j + bdiag[i+1]+1;
223b2b2dd24SShri Abhyankar     nz = bdiag[i] - bdiag[i+1] - 1;
224b2b2dd24SShri Abhyankar     for (j=0; j<nz; j++) {
2259566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2));
226b2b2dd24SShri Abhyankar     }
227b2b2dd24SShri Abhyankar   }
2289566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rtmp,mwork));
229ae3d28f0SHong Zhang 
2304dd39f65SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_2_NaturalOrdering;
2319f5c0bcdSBarry Smith   C->ops->forwardsolve   = MatForwardSolve_SeqBAIJ_2_NaturalOrdering;
2329f5c0bcdSBarry Smith   C->ops->backwardsolve  = MatBackwardSolve_SeqBAIJ_2_NaturalOrdering;
2334dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering;
234b2b2dd24SShri Abhyankar   C->assembled           = PETSC_TRUE;
23526fbe8dcSKarl Rupp 
2369566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333*2*2*2*n)); /* from inverting diagonal blocks */
237b2b2dd24SShri Abhyankar   PetscFunctionReturn(0);
238b2b2dd24SShri Abhyankar }
239b2b2dd24SShri Abhyankar 
24006e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_inplace(Mat B,Mat A,const MatFactorInfo *info)
2414eeb42bcSBarry Smith {
242719d5645SBarry Smith   Mat            C     = B;
2434eeb42bcSBarry Smith   Mat_SeqBAIJ    *a    = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
2447cdf2dbbSSatish Balay   IS             isrow = b->row,isicol = b->icol;
2455d0c19d7SBarry Smith   const PetscInt *r,*ic;
2465d0c19d7SBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
247690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row;
248690b6cddSBarry Smith   PetscInt       *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj;
249329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4;
2502515f8d2SSatish Balay   MatScalar      p1,p2,p3,p4;
2513f1db9ecSBarry Smith   MatScalar      *ba   = b->a,*aa = a->a;
252182b8fbaSHong Zhang   PetscReal      shift = info->shiftamount;
253a455e926SHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
2544eeb42bcSBarry Smith 
2553a40ed3dSBarry Smith   PetscFunctionBegin;
2560164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
2579566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow,&r));
2589566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol,&ic));
2599566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(4*(n+1),&rtmp));
2604eeb42bcSBarry Smith 
2614eeb42bcSBarry Smith   for (i=0; i<n; i++) {
2624078e994SBarry Smith     nz    = bi[i+1] - bi[i];
2634078e994SBarry Smith     ajtmp = bj + bi[i];
2644eeb42bcSBarry Smith     for  (j=0; j<nz; j++) {
2654eeb42bcSBarry Smith       x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0;
2664eeb42bcSBarry Smith     }
2674eeb42bcSBarry Smith     /* load in initial (unfactored row) */
2684eeb42bcSBarry Smith     idx      = r[i];
2694078e994SBarry Smith     nz       = ai[idx+1] - ai[idx];
2704078e994SBarry Smith     ajtmpold = aj + ai[idx];
2714078e994SBarry Smith     v        = aa + 4*ai[idx];
2724eeb42bcSBarry Smith     for (j=0; j<nz; j++) {
2734eeb42bcSBarry Smith       x    = rtmp+4*ic[ajtmpold[j]];
2744eeb42bcSBarry Smith       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
2754eeb42bcSBarry Smith       v   += 4;
2764eeb42bcSBarry Smith     }
2774eeb42bcSBarry Smith     row = *ajtmp++;
2784eeb42bcSBarry Smith     while (row < i) {
2794eeb42bcSBarry Smith       pc = rtmp + 4*row;
2804eeb42bcSBarry Smith       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
281d4a378daSJed Brown       if (p1 != (PetscScalar)0.0 || p2 != (PetscScalar)0.0 || p3 != (PetscScalar)0.0 || p4 != (PetscScalar)0.0) {
2824078e994SBarry Smith         pv    = ba + 4*diag_offset[row];
2834078e994SBarry Smith         pj    = bj + diag_offset[row] + 1;
2844eeb42bcSBarry Smith         x1    = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
2854eeb42bcSBarry Smith         pc[0] = m1 = p1*x1 + p3*x2;
2864eeb42bcSBarry Smith         pc[1] = m2 = p2*x1 + p4*x2;
2874eeb42bcSBarry Smith         pc[2] = m3 = p1*x3 + p3*x4;
2884eeb42bcSBarry Smith         pc[3] = m4 = p2*x3 + p4*x4;
2894078e994SBarry Smith         nz    = bi[row+1] - diag_offset[row] - 1;
2904eeb42bcSBarry Smith         pv   += 4;
2914eeb42bcSBarry Smith         for (j=0; j<nz; j++) {
2924eeb42bcSBarry Smith           x1    = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
2934eeb42bcSBarry Smith           x     = rtmp + 4*pj[j];
2944eeb42bcSBarry Smith           x[0] -= m1*x1 + m3*x2;
2954eeb42bcSBarry Smith           x[1] -= m2*x1 + m4*x2;
2964eeb42bcSBarry Smith           x[2] -= m1*x3 + m3*x4;
2974eeb42bcSBarry Smith           x[3] -= m2*x3 + m4*x4;
2984eeb42bcSBarry Smith           pv   += 4;
2994eeb42bcSBarry Smith         }
3009566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(16.0*nz+12.0));
3014eeb42bcSBarry Smith       }
3024eeb42bcSBarry Smith       row = *ajtmp++;
3034eeb42bcSBarry Smith     }
3044eeb42bcSBarry Smith     /* finished row so stick it into b->a */
3054078e994SBarry Smith     pv = ba + 4*bi[i];
3064078e994SBarry Smith     pj = bj + bi[i];
3074078e994SBarry Smith     nz = bi[i+1] - bi[i];
3084eeb42bcSBarry Smith     for (j=0; j<nz; j++) {
3094eeb42bcSBarry Smith       x     = rtmp+4*pj[j];
3104eeb42bcSBarry Smith       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
3114eeb42bcSBarry Smith       pv   += 4;
3124eeb42bcSBarry Smith     }
3134eeb42bcSBarry Smith     /* invert diagonal block */
3144078e994SBarry Smith     w    = ba + 4*diag_offset[i];
3159566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_2(w,shift,allowzeropivot,&zeropivotdetected));
3167b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3174eeb42bcSBarry Smith   }
3184eeb42bcSBarry Smith 
3199566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
3209566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol,&ic));
3219566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow,&r));
32226fbe8dcSKarl Rupp 
32306e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_2_inplace;
32406e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_inplace;
3254eeb42bcSBarry Smith   C->assembled           = PETSC_TRUE;
32626fbe8dcSKarl Rupp 
3279566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333*8*b->mbs)); /* from inverting diagonal blocks */
3283a40ed3dSBarry Smith   PetscFunctionReturn(0);
3294eeb42bcSBarry Smith }
3304cd38560SBarry Smith /*
3314cd38560SBarry Smith       Version for when blocks are 2 by 2 Using natural ordering
3324cd38560SBarry Smith */
33306e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
3344cd38560SBarry Smith {
3354cd38560SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
336690b6cddSBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
337690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row;
338690b6cddSBarry Smith   PetscInt       *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
339329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
3402515f8d2SSatish Balay   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4;
3414cd38560SBarry Smith   MatScalar      *ba   = b->a,*aa = a->a;
342182b8fbaSHong Zhang   PetscReal      shift = info->shiftamount;
343a455e926SHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
3444cd38560SBarry Smith 
3454cd38560SBarry Smith   PetscFunctionBegin;
3460164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
3479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(4*(n+1),&rtmp));
3484cd38560SBarry Smith   for (i=0; i<n; i++) {
3494cd38560SBarry Smith     nz    = bi[i+1] - bi[i];
3504cd38560SBarry Smith     ajtmp = bj + bi[i];
3514cd38560SBarry Smith     for  (j=0; j<nz; j++) {
3524cd38560SBarry Smith       x    = rtmp+4*ajtmp[j];
3534cd38560SBarry Smith       x[0] = x[1]  = x[2]  = x[3]  = 0.0;
3544cd38560SBarry Smith     }
3554cd38560SBarry Smith     /* load in initial (unfactored row) */
3564cd38560SBarry Smith     nz       = ai[i+1] - ai[i];
3574cd38560SBarry Smith     ajtmpold = aj + ai[i];
3584cd38560SBarry Smith     v        = aa + 4*ai[i];
3594cd38560SBarry Smith     for (j=0; j<nz; j++) {
3604cd38560SBarry Smith       x    = rtmp+4*ajtmpold[j];
3614cd38560SBarry Smith       x[0] = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
3624cd38560SBarry Smith       v   += 4;
3634cd38560SBarry Smith     }
3644cd38560SBarry Smith     row = *ajtmp++;
3654cd38560SBarry Smith     while (row < i) {
3664cd38560SBarry Smith       pc = rtmp + 4*row;
3674cd38560SBarry Smith       p1 = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
368d4a378daSJed Brown       if (p1 != (PetscScalar)0.0 || p2 != (PetscScalar)0.0 || p3 != (PetscScalar)0.0 || p4 != (PetscScalar)0.0) {
3694cd38560SBarry Smith         pv    = ba + 4*diag_offset[row];
3704cd38560SBarry Smith         pj    = bj + diag_offset[row] + 1;
3714cd38560SBarry Smith         x1    = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
3724cd38560SBarry Smith         pc[0] = m1 = p1*x1 + p3*x2;
3734cd38560SBarry Smith         pc[1] = m2 = p2*x1 + p4*x2;
3744cd38560SBarry Smith         pc[2] = m3 = p1*x3 + p3*x4;
3754cd38560SBarry Smith         pc[3] = m4 = p2*x3 + p4*x4;
3764cd38560SBarry Smith         nz    = bi[row+1] - diag_offset[row] - 1;
3774cd38560SBarry Smith         pv   += 4;
3784cd38560SBarry Smith         for (j=0; j<nz; j++) {
3794cd38560SBarry Smith           x1    = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
3804cd38560SBarry Smith           x     = rtmp + 4*pj[j];
3814cd38560SBarry Smith           x[0] -= m1*x1 + m3*x2;
3824cd38560SBarry Smith           x[1] -= m2*x1 + m4*x2;
3834cd38560SBarry Smith           x[2] -= m1*x3 + m3*x4;
3844cd38560SBarry Smith           x[3] -= m2*x3 + m4*x4;
3854cd38560SBarry Smith           pv   += 4;
3864cd38560SBarry Smith         }
3879566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(16.0*nz+12.0));
3884cd38560SBarry Smith       }
3894cd38560SBarry Smith       row = *ajtmp++;
3904cd38560SBarry Smith     }
3914cd38560SBarry Smith     /* finished row so stick it into b->a */
3924cd38560SBarry Smith     pv = ba + 4*bi[i];
3934cd38560SBarry Smith     pj = bj + bi[i];
3944cd38560SBarry Smith     nz = bi[i+1] - bi[i];
3954cd38560SBarry Smith     for (j=0; j<nz; j++) {
3964cd38560SBarry Smith       x     = rtmp+4*pj[j];
3974cd38560SBarry Smith       pv[0] = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
3982efa7f71SHong Zhang       /*
3992efa7f71SHong Zhang       printf(" col %d:",pj[j]);
4002efa7f71SHong Zhang       PetscInt j1;
4012efa7f71SHong Zhang       for (j1=0; j1<4; j1++) printf(" %g,",*(pv+j1));
4022efa7f71SHong Zhang       printf("\n");
4032efa7f71SHong Zhang       */
4044cd38560SBarry Smith       pv += 4;
4054cd38560SBarry Smith     }
4064cd38560SBarry Smith     /* invert diagonal block */
4074cd38560SBarry Smith     w = ba + 4*diag_offset[i];
4089566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_2(w,shift, allowzeropivot,&zeropivotdetected));
4097b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
4104cd38560SBarry Smith   }
4114cd38560SBarry Smith 
4129566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
41326fbe8dcSKarl Rupp 
41406e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_2_NaturalOrdering_inplace;
41506e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering_inplace;
4164cd38560SBarry Smith   C->assembled           = PETSC_TRUE;
41726fbe8dcSKarl Rupp 
4189566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333*8*b->mbs)); /* from inverting diagonal blocks */
4194cd38560SBarry Smith   PetscFunctionReturn(0);
4204cd38560SBarry Smith }
4217fc0212eSBarry Smith 
4227fc0212eSBarry Smith /* ----------------------------------------------------------- */
4237fc0212eSBarry Smith /*
4247fc0212eSBarry Smith      Version for when blocks are 1 by 1.
4257fc0212eSBarry Smith */
426048b5e81SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat B,Mat A,const MatFactorInfo *info)
427048b5e81SShri Abhyankar {
428048b5e81SShri Abhyankar   Mat             C     =B;
429048b5e81SShri Abhyankar   Mat_SeqBAIJ     *a    =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
430048b5e81SShri Abhyankar   IS              isrow = b->row,isicol = b->icol;
431048b5e81SShri Abhyankar   const PetscInt  *r,*ic,*ics;
432048b5e81SShri Abhyankar   const PetscInt  n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bdiag=b->diag;
433048b5e81SShri Abhyankar   PetscInt        i,j,k,nz,nzL,row,*pj;
434048b5e81SShri Abhyankar   const PetscInt  *ajtmp,*bjtmp;
435048b5e81SShri Abhyankar   MatScalar       *rtmp,*pc,multiplier,*pv;
436048b5e81SShri Abhyankar   const MatScalar *aa=a->a,*v;
437ace3abfcSBarry Smith   PetscBool       row_identity,col_identity;
438048b5e81SShri Abhyankar   FactorShiftCtx  sctx;
439048b5e81SShri Abhyankar   const PetscInt  *ddiag;
440048b5e81SShri Abhyankar   PetscReal       rs;
441048b5e81SShri Abhyankar   MatScalar       d;
442048b5e81SShri Abhyankar 
443048b5e81SShri Abhyankar   PetscFunctionBegin;
444048b5e81SShri Abhyankar   /* MatPivotSetUp(): initialize shift context sctx */
4459566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&sctx,sizeof(FactorShiftCtx)));
446048b5e81SShri Abhyankar 
447048b5e81SShri Abhyankar   if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
448048b5e81SShri Abhyankar     ddiag          = a->diag;
449048b5e81SShri Abhyankar     sctx.shift_top = info->zeropivot;
450048b5e81SShri Abhyankar     for (i=0; i<n; i++) {
451048b5e81SShri Abhyankar       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
452048b5e81SShri Abhyankar       d  = (aa)[ddiag[i]];
453048b5e81SShri Abhyankar       rs = -PetscAbsScalar(d) - PetscRealPart(d);
454048b5e81SShri Abhyankar       v  = aa+ai[i];
455048b5e81SShri Abhyankar       nz = ai[i+1] - ai[i];
45626fbe8dcSKarl Rupp       for (j=0; j<nz; j++) rs += PetscAbsScalar(v[j]);
457048b5e81SShri Abhyankar       if (rs>sctx.shift_top) sctx.shift_top = rs;
458048b5e81SShri Abhyankar     }
459048b5e81SShri Abhyankar     sctx.shift_top *= 1.1;
460048b5e81SShri Abhyankar     sctx.nshift_max = 5;
461048b5e81SShri Abhyankar     sctx.shift_lo   = 0.;
462048b5e81SShri Abhyankar     sctx.shift_hi   = 1.;
463048b5e81SShri Abhyankar   }
464048b5e81SShri Abhyankar 
4659566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow,&r));
4669566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol,&ic));
4679566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n+1,&rtmp));
468048b5e81SShri Abhyankar   ics  = ic;
469048b5e81SShri Abhyankar 
470048b5e81SShri Abhyankar   do {
471048b5e81SShri Abhyankar     sctx.newshift = PETSC_FALSE;
472048b5e81SShri Abhyankar     for (i=0; i<n; i++) {
473048b5e81SShri Abhyankar       /* zero rtmp */
474048b5e81SShri Abhyankar       /* L part */
475048b5e81SShri Abhyankar       nz    = bi[i+1] - bi[i];
476048b5e81SShri Abhyankar       bjtmp = bj + bi[i];
477048b5e81SShri Abhyankar       for  (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
478048b5e81SShri Abhyankar 
479048b5e81SShri Abhyankar       /* U part */
480048b5e81SShri Abhyankar       nz    = bdiag[i]-bdiag[i+1];
481048b5e81SShri Abhyankar       bjtmp = bj + bdiag[i+1]+1;
482048b5e81SShri Abhyankar       for  (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
483048b5e81SShri Abhyankar 
484048b5e81SShri Abhyankar       /* load in initial (unfactored row) */
485048b5e81SShri Abhyankar       nz    = ai[r[i]+1] - ai[r[i]];
486048b5e81SShri Abhyankar       ajtmp = aj + ai[r[i]];
487048b5e81SShri Abhyankar       v     = aa + ai[r[i]];
48826fbe8dcSKarl Rupp       for (j=0; j<nz; j++) rtmp[ics[ajtmp[j]]] = v[j];
48926fbe8dcSKarl Rupp 
490048b5e81SShri Abhyankar       /* ZeropivotApply() */
491048b5e81SShri Abhyankar       rtmp[i] += sctx.shift_amount;  /* shift the diagonal of the matrix */
492048b5e81SShri Abhyankar 
493048b5e81SShri Abhyankar       /* elimination */
494048b5e81SShri Abhyankar       bjtmp = bj + bi[i];
495048b5e81SShri Abhyankar       row   = *bjtmp++;
496048b5e81SShri Abhyankar       nzL   = bi[i+1] - bi[i];
497048b5e81SShri Abhyankar       for (k=0; k < nzL; k++) {
498048b5e81SShri Abhyankar         pc = rtmp + row;
499d4a378daSJed Brown         if (*pc != (PetscScalar)0.0) {
500048b5e81SShri Abhyankar           pv         = b->a + bdiag[row];
501048b5e81SShri Abhyankar           multiplier = *pc * (*pv);
502048b5e81SShri Abhyankar           *pc        = multiplier;
50326fbe8dcSKarl Rupp 
504048b5e81SShri Abhyankar           pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
505048b5e81SShri Abhyankar           pv = b->a + bdiag[row+1]+1;
506048b5e81SShri Abhyankar           nz = bdiag[row]-bdiag[row+1]-1; /* num of entries in U(row,:) excluding diag */
507048b5e81SShri Abhyankar           for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
5089566063dSJacob Faibussowitsch           PetscCall(PetscLogFlops(2.0*nz));
509048b5e81SShri Abhyankar         }
510048b5e81SShri Abhyankar         row = *bjtmp++;
511048b5e81SShri Abhyankar       }
512048b5e81SShri Abhyankar 
513048b5e81SShri Abhyankar       /* finished row so stick it into b->a */
514048b5e81SShri Abhyankar       rs = 0.0;
515048b5e81SShri Abhyankar       /* L part */
516048b5e81SShri Abhyankar       pv = b->a + bi[i];
517048b5e81SShri Abhyankar       pj = b->j + bi[i];
518048b5e81SShri Abhyankar       nz = bi[i+1] - bi[i];
519048b5e81SShri Abhyankar       for (j=0; j<nz; j++) {
520048b5e81SShri Abhyankar         pv[j] = rtmp[pj[j]]; rs += PetscAbsScalar(pv[j]);
521048b5e81SShri Abhyankar       }
522048b5e81SShri Abhyankar 
523048b5e81SShri Abhyankar       /* U part */
524048b5e81SShri Abhyankar       pv = b->a + bdiag[i+1]+1;
525048b5e81SShri Abhyankar       pj = b->j + bdiag[i+1]+1;
526048b5e81SShri Abhyankar       nz = bdiag[i] - bdiag[i+1]-1;
527048b5e81SShri Abhyankar       for (j=0; j<nz; j++) {
528048b5e81SShri Abhyankar         pv[j] = rtmp[pj[j]]; rs += PetscAbsScalar(pv[j]);
529048b5e81SShri Abhyankar       }
530048b5e81SShri Abhyankar 
531048b5e81SShri Abhyankar       sctx.rs = rs;
532048b5e81SShri Abhyankar       sctx.pv = rtmp[i];
5339566063dSJacob Faibussowitsch       PetscCall(MatPivotCheck(B,A,info,&sctx,i));
534048b5e81SShri Abhyankar       if (sctx.newshift) break; /* break for-loop */
535048b5e81SShri Abhyankar       rtmp[i] = sctx.pv; /* sctx.pv might be updated in the case of MAT_SHIFT_INBLOCKS */
536048b5e81SShri Abhyankar 
537a5b23f4aSJose E. Roman       /* Mark diagonal and invert diagonal for simpler triangular solves */
538048b5e81SShri Abhyankar       pv  = b->a + bdiag[i];
539d4a378daSJed Brown       *pv = (PetscScalar)1.0/rtmp[i];
540048b5e81SShri Abhyankar 
541048b5e81SShri Abhyankar     } /* endof for (i=0; i<n; i++) { */
542048b5e81SShri Abhyankar 
543048b5e81SShri Abhyankar     /* MatPivotRefine() */
544048b5e81SShri Abhyankar     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE && !sctx.newshift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
545048b5e81SShri Abhyankar       /*
546048b5e81SShri Abhyankar        * if no shift in this attempt & shifting & started shifting & can refine,
547048b5e81SShri Abhyankar        * then try lower shift
548048b5e81SShri Abhyankar        */
549048b5e81SShri Abhyankar       sctx.shift_hi       = sctx.shift_fraction;
550048b5e81SShri Abhyankar       sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
551048b5e81SShri Abhyankar       sctx.shift_amount   = sctx.shift_fraction * sctx.shift_top;
552048b5e81SShri Abhyankar       sctx.newshift       = PETSC_TRUE;
553048b5e81SShri Abhyankar       sctx.nshift++;
554048b5e81SShri Abhyankar     }
555048b5e81SShri Abhyankar   } while (sctx.newshift);
556048b5e81SShri Abhyankar 
5579566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
5589566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol,&ic));
5599566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow,&r));
560048b5e81SShri Abhyankar 
5619566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isrow,&row_identity));
5629566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isicol,&col_identity));
563048b5e81SShri Abhyankar   if (row_identity && col_identity) {
564048b5e81SShri Abhyankar     C->ops->solve          = MatSolve_SeqBAIJ_1_NaturalOrdering;
5659f5c0bcdSBarry Smith     C->ops->forwardsolve   = MatForwardSolve_SeqBAIJ_1_NaturalOrdering;
5669f5c0bcdSBarry Smith     C->ops->backwardsolve  = MatBackwardSolve_SeqBAIJ_1_NaturalOrdering;
56793fd935bSShri Abhyankar     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering;
568048b5e81SShri Abhyankar   } else {
569048b5e81SShri Abhyankar     C->ops->solve          = MatSolve_SeqBAIJ_1;
57093fd935bSShri Abhyankar     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1;
571048b5e81SShri Abhyankar   }
572048b5e81SShri Abhyankar   C->assembled = PETSC_TRUE;
5739566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(C->cmap->n));
574048b5e81SShri Abhyankar 
575048b5e81SShri Abhyankar   /* MatShiftView(A,info,&sctx) */
576048b5e81SShri Abhyankar   if (sctx.nshift) {
577048b5e81SShri Abhyankar     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
5789566063dSJacob 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));
579048b5e81SShri Abhyankar     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
5809566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A,"number of shift_nz tries %" PetscInt_FMT ", shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount));
581048b5e81SShri Abhyankar     } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) {
5829566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A,"number of shift_inblocks applied %" PetscInt_FMT ", each shift_amount %g\n",sctx.nshift,(double)info->shiftamount));
583048b5e81SShri Abhyankar     }
584048b5e81SShri Abhyankar   }
585048b5e81SShri Abhyankar   PetscFunctionReturn(0);
586048b5e81SShri Abhyankar }
587048b5e81SShri Abhyankar 
58806e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info)
5897fc0212eSBarry Smith {
5907fc0212eSBarry Smith   Mat_SeqBAIJ    *a    = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
5917cdf2dbbSSatish Balay   IS             isrow = b->row,isicol = b->icol;
5925d0c19d7SBarry Smith   const PetscInt *r,*ic;
5935d0c19d7SBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
594690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j;
595690b6cddSBarry Smith   PetscInt       *diag_offset = b->diag,diag,*pj;
596329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,multiplier,*pc;
5973f1db9ecSBarry Smith   MatScalar      *ba = b->a,*aa = a->a;
598ace3abfcSBarry Smith   PetscBool      row_identity, col_identity;
5997fc0212eSBarry Smith 
6003a40ed3dSBarry Smith   PetscFunctionBegin;
6019566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow,&r));
6029566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol,&ic));
6039566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n+1,&rtmp));
6047fc0212eSBarry Smith 
6057fc0212eSBarry Smith   for (i=0; i<n; i++) {
6064078e994SBarry Smith     nz    = bi[i+1] - bi[i];
6074078e994SBarry Smith     ajtmp = bj + bi[i];
6087fc0212eSBarry Smith     for  (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0;
6097fc0212eSBarry Smith 
6107fc0212eSBarry Smith     /* load in initial (unfactored row) */
6114078e994SBarry Smith     nz       = ai[r[i]+1] - ai[r[i]];
6124078e994SBarry Smith     ajtmpold = aj + ai[r[i]];
6134078e994SBarry Smith     v        = aa + ai[r[i]];
6147fc0212eSBarry Smith     for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] =  v[j];
6157fc0212eSBarry Smith 
6167fc0212eSBarry Smith     row = *ajtmp++;
6177fc0212eSBarry Smith     while (row < i) {
6187fc0212eSBarry Smith       pc = rtmp + row;
6197fc0212eSBarry Smith       if (*pc != 0.0) {
6204078e994SBarry Smith         pv         = ba + diag_offset[row];
6214078e994SBarry Smith         pj         = bj + diag_offset[row] + 1;
6227fc0212eSBarry Smith         multiplier = *pc * *pv++;
6237fc0212eSBarry Smith         *pc        = multiplier;
6244078e994SBarry Smith         nz         = bi[row+1] - diag_offset[row] - 1;
6257fc0212eSBarry Smith         for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
6269566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(1.0+2.0*nz));
6277fc0212eSBarry Smith       }
6287fc0212eSBarry Smith       row = *ajtmp++;
6297fc0212eSBarry Smith     }
6307fc0212eSBarry Smith     /* finished row so stick it into b->a */
6314078e994SBarry Smith     pv = ba + bi[i];
6324078e994SBarry Smith     pj = bj + bi[i];
6334078e994SBarry Smith     nz = bi[i+1] - bi[i];
63426fbe8dcSKarl Rupp     for (j=0; j<nz; j++) pv[j] = rtmp[pj[j]];
6354078e994SBarry Smith     diag = diag_offset[i] - bi[i];
6367fc0212eSBarry Smith     /* check pivot entry for current row */
63708401ef6SPierre 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);
6387fc0212eSBarry Smith     pv[diag] = 1.0/pv[diag];
6397fc0212eSBarry Smith   }
6407fc0212eSBarry Smith 
6419566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
6429566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol,&ic));
6439566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow,&r));
6449566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isrow,&row_identity));
6459566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isicol,&col_identity));
646f58c8c74SBarry Smith   if (row_identity && col_identity) {
64706e38f1dSHong Zhang     C->ops->solve          = MatSolve_SeqBAIJ_1_NaturalOrdering_inplace;
64806e38f1dSHong Zhang     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering_inplace;
649f58c8c74SBarry Smith   } else {
65006e38f1dSHong Zhang     C->ops->solve          = MatSolve_SeqBAIJ_1_inplace;
65106e38f1dSHong Zhang     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_inplace;
652f58c8c74SBarry Smith   }
6537fc0212eSBarry Smith   C->assembled = PETSC_TRUE;
6549566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(C->cmap->n));
6553a40ed3dSBarry Smith   PetscFunctionReturn(0);
6567fc0212eSBarry Smith }
6577fc0212eSBarry Smith 
6584ac6704cSBarry Smith static PetscErrorCode MatFactorGetSolverType_petsc(Mat A,MatSolverType *type)
6594ac6704cSBarry Smith {
6604ac6704cSBarry Smith   PetscFunctionBegin;
6614ac6704cSBarry Smith   *type = MATSOLVERPETSC;
6624ac6704cSBarry Smith   PetscFunctionReturn(0);
6634ac6704cSBarry Smith }
6644ac6704cSBarry Smith 
665fd9d3c67SJed Brown PETSC_INTERN PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A,MatFactorType ftype,Mat *B)
666b24902e0SBarry Smith {
667d0f46423SBarry Smith   PetscInt       n = A->rmap->n;
668b24902e0SBarry Smith 
669b24902e0SBarry Smith   PetscFunctionBegin;
670b499a2aeSHong Zhang #if defined(PETSC_USE_COMPLEX)
671*b94d7dedSBarry 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");
672b499a2aeSHong Zhang #endif
6739566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A),B));
6749566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*B,n,n,n,n));
675c8342467SHong Zhang   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) {
6769566063dSJacob Faibussowitsch     PetscCall(MatSetType(*B,MATSEQBAIJ));
67726fbe8dcSKarl Rupp 
6784dd39f65SShri Abhyankar     (*B)->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqBAIJ;
6794dd39f65SShri Abhyankar     (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqBAIJ;
6809566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGND,(char**)&(*B)->preferredordering[MAT_FACTOR_LU]));
6819566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGNATURAL,(char**)&(*B)->preferredordering[MAT_FACTOR_ILU]));
6829566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGNATURAL,(char**)&(*B)->preferredordering[MAT_FACTOR_ILUDT]));
683b24902e0SBarry Smith   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
6849566063dSJacob Faibussowitsch     PetscCall(MatSetType(*B,MATSEQSBAIJ));
6859566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(*B,A->rmap->bs,MAT_SKIP_ALLOCATION,NULL));
68626fbe8dcSKarl Rupp 
6875c9eb25fSBarry Smith     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqBAIJ;
6885c9eb25fSBarry Smith     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ;
6894ac6704cSBarry Smith     /*  Future optimization would be direct symbolic and numerical factorization for BAIJ to support orderings and Cholesky, instead of first converting to SBAIJ */
6909566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGNATURAL,(char**)&(*B)->preferredordering[MAT_FACTOR_CHOLESKY]));
6919566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGNATURAL,(char**)&(*B)->preferredordering[MAT_FACTOR_ICC]));
692e32f2f54SBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
693d5f3da31SBarry Smith   (*B)->factortype = ftype;
694f73b0415SBarry Smith   (*B)->canuseordering = PETSC_TRUE;
69500c67f3bSHong Zhang 
6969566063dSJacob Faibussowitsch   PetscCall(PetscFree((*B)->solvertype));
6979566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERPETSC,&(*B)->solvertype));
6989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)*B,"MatFactorGetSolverType_C",MatFactorGetSolverType_petsc));
699b24902e0SBarry Smith   PetscFunctionReturn(0);
700b24902e0SBarry Smith }
701273d9f13SBarry Smith 
7025d517e7eSBarry Smith /* ----------------------------------------------------------- */
7030481f469SBarry Smith PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,const MatFactorInfo *info)
7045d517e7eSBarry Smith {
7055d517e7eSBarry Smith   Mat            C;
7065d517e7eSBarry Smith 
7073a40ed3dSBarry Smith   PetscFunctionBegin;
7089566063dSJacob Faibussowitsch   PetscCall(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&C));
7099566063dSJacob Faibussowitsch   PetscCall(MatLUFactorSymbolic(C,A,row,col,info));
7109566063dSJacob Faibussowitsch   PetscCall(MatLUFactorNumeric(C,A,info));
71126fbe8dcSKarl Rupp 
712db4efbfdSBarry Smith   A->ops->solve          = C->ops->solve;
713db4efbfdSBarry Smith   A->ops->solvetranspose = C->ops->solvetranspose;
71426fbe8dcSKarl Rupp 
7159566063dSJacob Faibussowitsch   PetscCall(MatHeaderMerge(A,&C));
7169566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)A,(PetscObject)((Mat_SeqBAIJ*)(A->data))->icol));
7173a40ed3dSBarry Smith   PetscFunctionReturn(0);
7185d517e7eSBarry Smith }
7194cd38560SBarry Smith 
720c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h>
7210481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
722c05c3958SHong Zhang {
72378910aadSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
72478910aadSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
72578910aadSHong Zhang   IS             ip=b->row;
7265d0c19d7SBarry Smith   const PetscInt *rip;
7275d0c19d7SBarry Smith   PetscInt       i,j,mbs=a->mbs,bs=A->rmap->bs,*bi=b->i,*bj=b->j,*bcol;
72878910aadSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
72978910aadSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
73078910aadSHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
73107b50cabSHong Zhang   PetscReal      rs;
7320e95ead3SHong Zhang   FactorShiftCtx sctx;
73378910aadSHong Zhang 
734c05c3958SHong Zhang   PetscFunctionBegin;
73507b50cabSHong Zhang   if (bs > 1) { /* convert A to a SBAIJ matrix and apply Cholesky factorization from it */
7366ad2eaddSHong Zhang     if (!a->sbaijMat) {
7379566063dSJacob Faibussowitsch       PetscCall(MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat));
7386ad2eaddSHong Zhang     }
7399566063dSJacob Faibussowitsch     PetscCall((a->sbaijMat)->ops->choleskyfactornumeric(C,a->sbaijMat,info));
7409566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&a->sbaijMat));
7416ad2eaddSHong Zhang     PetscFunctionReturn(0);
7426ad2eaddSHong Zhang   }
74378910aadSHong Zhang 
74407b50cabSHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
7459566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&sctx,sizeof(FactorShiftCtx)));
7463cea4cbeSHong Zhang 
7479566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(ip,&rip));
7489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&jl));
74978910aadSHong Zhang 
75075567043SBarry Smith   sctx.shift_amount = 0.;
7513cea4cbeSHong Zhang   sctx.nshift       = 0;
75278910aadSHong Zhang   do {
75307b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
75478910aadSHong Zhang     for (i=0; i<mbs; i++) {
755e60cf9a0SBarry Smith       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
75678910aadSHong Zhang     }
75778910aadSHong Zhang 
75878910aadSHong Zhang     for (k = 0; k<mbs; k++) {
75978910aadSHong Zhang       bval = ba + bi[k];
76078910aadSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
76178910aadSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
76278910aadSHong Zhang       for (j = jmin; j < jmax; j++) {
76378910aadSHong Zhang         col = rip[aj[j]];
76478910aadSHong Zhang         if (col >= k) { /* only take upper triangular entry */
76578910aadSHong Zhang           rtmp[col] = aa[j];
76678910aadSHong Zhang           *bval++   = 0.0; /* for in-place factorization */
76778910aadSHong Zhang         }
76878910aadSHong Zhang       }
7693cea4cbeSHong Zhang 
7703cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
7713cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
77278910aadSHong Zhang 
77378910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
77478910aadSHong Zhang       dk = rtmp[k];
77578910aadSHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
77678910aadSHong Zhang 
77778910aadSHong Zhang       while (i < k) {
77878910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
77978910aadSHong Zhang 
78078910aadSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
78178910aadSHong Zhang         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
78278910aadSHong Zhang         uikdi   = -ba[ili]*ba[bi[i]]; /* diagonal(k) */
78378910aadSHong Zhang         dk     += uikdi*ba[ili];
78478910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
78578910aadSHong Zhang 
78678910aadSHong Zhang         /* add multiple of row i to k-th row */
78778910aadSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
78878910aadSHong Zhang         if (jmin < jmax) {
78978910aadSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
79078910aadSHong Zhang           /* update il and jl for row i */
79178910aadSHong Zhang           il[i] = jmin;
79278910aadSHong Zhang           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
79378910aadSHong Zhang         }
79478910aadSHong Zhang         i = nexti;
79578910aadSHong Zhang       }
79678910aadSHong Zhang 
7973cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
7983cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
7993cea4cbeSHong Zhang       rs   = 0.0;
80078910aadSHong Zhang       jmin = bi[k]+1;
80178910aadSHong Zhang       nz   = bi[k+1] - jmin;
80278910aadSHong Zhang       if (nz) {
80378910aadSHong Zhang         bcol = bj + jmin;
80478910aadSHong Zhang         while (nz--) {
80589f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
80689f845c8SHong Zhang           bcol++;
80778910aadSHong Zhang         }
80878910aadSHong Zhang       }
8093cea4cbeSHong Zhang 
8103cea4cbeSHong Zhang       sctx.rs = rs;
8113cea4cbeSHong Zhang       sctx.pv = dk;
8129566063dSJacob Faibussowitsch       PetscCall(MatPivotCheck(C,A,info,&sctx,k));
81307b50cabSHong Zhang       if (sctx.newshift) break;
8140e95ead3SHong Zhang       dk = sctx.pv;
81578910aadSHong Zhang 
81678910aadSHong Zhang       /* copy data into U(k,:) */
81778910aadSHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
81878910aadSHong Zhang       jmin      = bi[k]+1; jmax = bi[k+1];
81978910aadSHong Zhang       if (jmin < jmax) {
82078910aadSHong Zhang         for (j=jmin; j<jmax; j++) {
82178910aadSHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
82278910aadSHong Zhang         }
82378910aadSHong Zhang         /* add the k-th row into il and jl */
82478910aadSHong Zhang         il[k] = jmin;
82578910aadSHong Zhang         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
82678910aadSHong Zhang       }
82778910aadSHong Zhang     }
82807b50cabSHong Zhang   } while (sctx.newshift);
8299566063dSJacob Faibussowitsch   PetscCall(PetscFree3(rtmp,il,jl));
83078910aadSHong Zhang 
8319566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(ip,&rip));
83226fbe8dcSKarl Rupp 
83378910aadSHong Zhang   C->assembled    = PETSC_TRUE;
83478910aadSHong Zhang   C->preallocated = PETSC_TRUE;
83526fbe8dcSKarl Rupp 
8369566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(C->rmap->N));
8373cea4cbeSHong Zhang   if (sctx.nshift) {
838f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
8399566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A,"number of shiftpd tries %" PetscInt_FMT ", shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount));
840f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
8419566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A,"number of shiftnz tries %" PetscInt_FMT ", shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount));
84278910aadSHong Zhang     }
84378910aadSHong Zhang   }
844c05c3958SHong Zhang   PetscFunctionReturn(0);
845c05c3958SHong Zhang }
8464cd38560SBarry Smith 
8470481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
848c05c3958SHong Zhang {
84978910aadSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
85078910aadSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
8513cea4cbeSHong Zhang   PetscInt       i,j,am=a->mbs;
85278910aadSHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
8533cea4cbeSHong Zhang   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
85478910aadSHong Zhang   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
8550e95ead3SHong Zhang   PetscReal      rs;
8560e95ead3SHong Zhang   FactorShiftCtx sctx;
85778910aadSHong Zhang 
858c05c3958SHong Zhang   PetscFunctionBegin;
8590e95ead3SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
8609566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&sctx,sizeof(FactorShiftCtx)));
8613cea4cbeSHong Zhang 
8629566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(am,&rtmp,am,&il,am,&jl));
86378910aadSHong Zhang 
86478910aadSHong Zhang   do {
86507b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
86678910aadSHong Zhang     for (i=0; i<am; i++) {
867e60cf9a0SBarry Smith       rtmp[i] = 0.0; jl[i] = am; il[0] = 0;
86878910aadSHong Zhang     }
86978910aadSHong Zhang 
87078910aadSHong Zhang     for (k = 0; k<am; k++) {
87178910aadSHong Zhang       /* initialize k-th row with elements nonzero in row perm(k) of A */
87278910aadSHong Zhang       nz   = ai[k+1] - ai[k];
87378910aadSHong Zhang       acol = aj + ai[k];
87478910aadSHong Zhang       aval = aa + ai[k];
87578910aadSHong Zhang       bval = ba + bi[k];
87678910aadSHong Zhang       while (nz--) {
87778910aadSHong Zhang         if (*acol < k) { /* skip lower triangular entries */
87878910aadSHong Zhang           acol++; aval++;
87978910aadSHong Zhang         } else {
88078910aadSHong Zhang           rtmp[*acol++] = *aval++;
88178910aadSHong Zhang           *bval++       = 0.0; /* for in-place factorization */
88278910aadSHong Zhang         }
88378910aadSHong Zhang       }
8843cea4cbeSHong Zhang 
8853cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
8863cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
88778910aadSHong Zhang 
88878910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
88978910aadSHong Zhang       dk = rtmp[k];
89078910aadSHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
89178910aadSHong Zhang 
89278910aadSHong Zhang       while (i < k) {
89378910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
89478910aadSHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
89578910aadSHong Zhang         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
89678910aadSHong Zhang         uikdi   = -ba[ili]*ba[bi[i]];
89778910aadSHong Zhang         dk     += uikdi*ba[ili];
89878910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
89978910aadSHong Zhang 
90078910aadSHong Zhang         /* add multiple of row i to k-th row ... */
90178910aadSHong Zhang         jmin = ili + 1;
90278910aadSHong Zhang         nz   = bi[i+1] - jmin;
90378910aadSHong Zhang         if (nz > 0) {
90478910aadSHong Zhang           bcol = bj + jmin;
90578910aadSHong Zhang           bval = ba + jmin;
90678910aadSHong Zhang           while (nz--) rtmp[*bcol++] += uikdi*(*bval++);
90778910aadSHong Zhang           /* update il and jl for i-th row */
90878910aadSHong Zhang           il[i] = jmin;
90978910aadSHong Zhang           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
91078910aadSHong Zhang         }
91178910aadSHong Zhang         i = nexti;
91278910aadSHong Zhang       }
91378910aadSHong Zhang 
9143cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
9153cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
9163cea4cbeSHong Zhang       rs   = 0.0;
91778910aadSHong Zhang       jmin = bi[k]+1;
91878910aadSHong Zhang       nz   = bi[k+1] - jmin;
91978910aadSHong Zhang       if (nz) {
92078910aadSHong Zhang         bcol = bj + jmin;
92178910aadSHong Zhang         while (nz--) {
92289f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
92389f845c8SHong Zhang           bcol++;
92478910aadSHong Zhang         }
92578910aadSHong Zhang       }
9263cea4cbeSHong Zhang 
9273cea4cbeSHong Zhang       sctx.rs = rs;
9283cea4cbeSHong Zhang       sctx.pv = dk;
9299566063dSJacob Faibussowitsch       PetscCall(MatPivotCheck(C,A,info,&sctx,k));
93007b50cabSHong Zhang       if (sctx.newshift) break;    /* sctx.shift_amount is updated */
9310e95ead3SHong Zhang       dk = sctx.pv;
93278910aadSHong Zhang 
93378910aadSHong Zhang       /* copy data into U(k,:) */
93478910aadSHong Zhang       ba[bi[k]] = 1.0/dk;
93578910aadSHong Zhang       jmin      = bi[k]+1;
93678910aadSHong Zhang       nz        = bi[k+1] - jmin;
93778910aadSHong Zhang       if (nz) {
93878910aadSHong Zhang         bcol = bj + jmin;
93978910aadSHong Zhang         bval = ba + jmin;
94078910aadSHong Zhang         while (nz--) {
94178910aadSHong Zhang           *bval++       = rtmp[*bcol];
94278910aadSHong Zhang           rtmp[*bcol++] = 0.0;
94378910aadSHong Zhang         }
94478910aadSHong Zhang         /* add k-th row into il and jl */
94578910aadSHong Zhang         il[k] = jmin;
94678910aadSHong Zhang         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
94778910aadSHong Zhang       }
94878910aadSHong Zhang     }
94907b50cabSHong Zhang   } while (sctx.newshift);
9509566063dSJacob Faibussowitsch   PetscCall(PetscFree3(rtmp,il,jl));
95178910aadSHong Zhang 
9520a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
9530a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
95478910aadSHong Zhang   C->assembled           = PETSC_TRUE;
95578910aadSHong Zhang   C->preallocated        = PETSC_TRUE;
95626fbe8dcSKarl Rupp 
9579566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(C->rmap->N));
9583cea4cbeSHong Zhang   if (sctx.nshift) {
959f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
9609566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A,"number of shiftnz tries %" PetscInt_FMT ", shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount));
961f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
9629566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A,"number of shiftpd tries %" PetscInt_FMT ", shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount));
96378910aadSHong Zhang     }
96478910aadSHong Zhang   }
965c05c3958SHong Zhang   PetscFunctionReturn(0);
966c05c3958SHong Zhang }
967c05c3958SHong Zhang 
968c6db04a5SJed Brown #include <petscbt.h>
969c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
9700481f469SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
971c05c3958SHong Zhang {
97278910aadSHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
97378910aadSHong Zhang   Mat_SeqSBAIJ       *b;
97478910aadSHong Zhang   Mat                B;
97548dd3d27SHong Zhang   PetscBool          perm_identity,missing;
9765d0c19d7SBarry Smith   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap->bs,*ui;
9775d0c19d7SBarry Smith   const PetscInt     *rip;
97878910aadSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
9790298fd71SBarry Smith   PetscInt           nlnk,*lnk,*lnk_lvl=NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
98078910aadSHong Zhang   PetscReal          fill          =info->fill,levels=info->levels;
9810298fd71SBarry Smith   PetscFreeSpaceList free_space    =NULL,current_space=NULL;
9820298fd71SBarry Smith   PetscFreeSpaceList free_space_lvl=NULL,current_space_lvl=NULL;
98378910aadSHong Zhang   PetscBT            lnkbt;
98478910aadSHong Zhang 
985c05c3958SHong Zhang   PetscFunctionBegin;
9869566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal(A,&missing,&i));
98728b400f6SJacob Faibussowitsch   PetscCheck(!missing,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %" PetscInt_FMT,i);
98848dd3d27SHong Zhang 
9896ad2eaddSHong Zhang   if (bs > 1) {
9906ad2eaddSHong Zhang     if (!a->sbaijMat) {
9919566063dSJacob Faibussowitsch       PetscCall(MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat));
9926ad2eaddSHong Zhang     }
993719d5645SBarry Smith     (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ;  /* undue the change made in MatGetFactor_seqbaij_petsc */
99426fbe8dcSKarl Rupp 
9959566063dSJacob Faibussowitsch     PetscCall(MatICCFactorSymbolic(fact,a->sbaijMat,perm,info));
9966ad2eaddSHong Zhang     PetscFunctionReturn(0);
9976ad2eaddSHong Zhang   }
9986ad2eaddSHong Zhang 
9999566063dSJacob Faibussowitsch   PetscCall(ISIdentity(perm,&perm_identity));
10009566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(perm,&rip));
100178910aadSHong Zhang 
100278910aadSHong Zhang   /* special case that simply copies fill pattern */
100378910aadSHong Zhang   if (!levels && perm_identity) {
10049566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(am+1,&ui));
100526fbe8dcSKarl Rupp     for (i=0; i<am; i++) ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
1006719d5645SBarry Smith     B    = fact;
10079566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(B,1,0,ui));
100878910aadSHong Zhang 
100978910aadSHong Zhang     b  = (Mat_SeqSBAIJ*)B->data;
101078910aadSHong Zhang     uj = b->j;
101178910aadSHong Zhang     for (i=0; i<am; i++) {
101278910aadSHong Zhang       aj = a->j + a->diag[i];
101326fbe8dcSKarl Rupp       for (j=0; j<ui[i]; j++) *uj++ = *aj++;
101478910aadSHong Zhang       b->ilen[i] = ui[i];
101578910aadSHong Zhang     }
10169566063dSJacob Faibussowitsch     PetscCall(PetscFree(ui));
101726fbe8dcSKarl Rupp 
1018d5f3da31SBarry Smith     B->factortype = MAT_FACTOR_NONE;
101926fbe8dcSKarl Rupp 
10209566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
10219566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
1022d5f3da31SBarry Smith     B->factortype = MAT_FACTOR_ICC;
102378910aadSHong Zhang 
102478910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
102578910aadSHong Zhang     PetscFunctionReturn(0);
102678910aadSHong Zhang   }
102778910aadSHong Zhang 
102878910aadSHong Zhang   /* initialization */
10299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(am+1,&ui));
1030e60cf9a0SBarry Smith   ui[0] = 0;
10319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2*am+1,&cols_lvl));
103278910aadSHong Zhang 
103378910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
103478910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
10359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(am,&uj_ptr,am,&uj_lvl_ptr,am,&il,am,&jl));
103678910aadSHong Zhang   for (i=0; i<am; i++) {
1037e60cf9a0SBarry Smith     jl[i] = am; il[i] = 0;
103878910aadSHong Zhang   }
103978910aadSHong Zhang 
104078910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
104178910aadSHong Zhang   nlnk = am + 1;
10429566063dSJacob Faibussowitsch   PetscCall(PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt));
104378910aadSHong Zhang 
104495b5ac22SHong Zhang   /* initial FreeSpace size is fill*(ai[am]+am)/2 */
10459566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am]/2,am/2)),&free_space));
104626fbe8dcSKarl Rupp 
104778910aadSHong Zhang   current_space = free_space;
104826fbe8dcSKarl Rupp 
10499566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am]/2,am/2)),&free_space_lvl));
105078910aadSHong Zhang   current_space_lvl = free_space_lvl;
105178910aadSHong Zhang 
105278910aadSHong Zhang   for (k=0; k<am; k++) {  /* for each active row k */
105378910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
105478910aadSHong Zhang     nzk         = 0;
105578910aadSHong Zhang     ncols       = ai[rip[k]+1] - ai[rip[k]];
105678910aadSHong Zhang     ncols_upper = 0;
105778910aadSHong Zhang     cols        = cols_lvl + am;
105878910aadSHong Zhang     for (j=0; j<ncols; j++) {
105978910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
106078910aadSHong Zhang       if (i >= k) { /* only take upper triangular entry */
106178910aadSHong Zhang         cols[ncols_upper]     = i;
106278910aadSHong Zhang         cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
106378910aadSHong Zhang         ncols_upper++;
106478910aadSHong Zhang       }
106578910aadSHong Zhang     }
10669566063dSJacob Faibussowitsch     PetscCall(PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,&nlnk,lnk,lnk_lvl,lnkbt));
106778910aadSHong Zhang     nzk += nlnk;
106878910aadSHong Zhang 
106978910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
107078910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
107178910aadSHong Zhang 
107278910aadSHong Zhang     while (prow < k) {
107378910aadSHong Zhang       nextprow = jl[prow];
107478910aadSHong Zhang 
107578910aadSHong Zhang       /* merge prow into k-th row */
107678910aadSHong Zhang       jmin  = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */
107778910aadSHong Zhang       jmax  = ui[prow+1];
107878910aadSHong Zhang       ncols = jmax-jmin;
107978910aadSHong Zhang       i     = jmin - ui[prow];
108078910aadSHong Zhang       cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
108178910aadSHong Zhang       for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
10829566063dSJacob Faibussowitsch       PetscCall(PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,&nlnk,lnk,lnk_lvl,lnkbt));
108378910aadSHong Zhang       nzk += nlnk;
108478910aadSHong Zhang 
108578910aadSHong Zhang       /* update il and jl for prow */
108678910aadSHong Zhang       if (jmin < jmax) {
108778910aadSHong Zhang         il[prow] = jmin;
108826fbe8dcSKarl Rupp 
108978910aadSHong Zhang         j = *cols; jl[prow] = jl[j]; jl[j] = prow;
109078910aadSHong Zhang       }
109178910aadSHong Zhang       prow = nextprow;
109278910aadSHong Zhang     }
109378910aadSHong Zhang 
109478910aadSHong Zhang     /* if free space is not available, make more free space */
109578910aadSHong Zhang     if (current_space->local_remaining<nzk) {
109678910aadSHong Zhang       i    = am - k + 1; /* num of unfactored rows */
1097f91af8c7SBarry Smith       i    = PetscMin(PetscIntMultTruncate(i,nzk), PetscIntMultTruncate(i,i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
10989566063dSJacob Faibussowitsch       PetscCall(PetscFreeSpaceGet(i,&current_space));
10999566063dSJacob Faibussowitsch       PetscCall(PetscFreeSpaceGet(i,&current_space_lvl));
110078910aadSHong Zhang       reallocs++;
110178910aadSHong Zhang     }
110278910aadSHong Zhang 
110378910aadSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
11049566063dSJacob Faibussowitsch     PetscCall(PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt));
110578910aadSHong Zhang 
110678910aadSHong Zhang     /* add the k-th row into il and jl */
110778910aadSHong Zhang     if (nzk-1 > 0) {
110878910aadSHong Zhang       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
110978910aadSHong Zhang       jl[k] = jl[i]; jl[i] = k;
111078910aadSHong Zhang       il[k] = ui[k] + 1;
111178910aadSHong Zhang     }
111278910aadSHong Zhang     uj_ptr[k]     = current_space->array;
111378910aadSHong Zhang     uj_lvl_ptr[k] = current_space_lvl->array;
111478910aadSHong Zhang 
111578910aadSHong Zhang     current_space->array           += nzk;
111678910aadSHong Zhang     current_space->local_used      += nzk;
111778910aadSHong Zhang     current_space->local_remaining -= nzk;
111878910aadSHong Zhang 
111978910aadSHong Zhang     current_space_lvl->array           += nzk;
112078910aadSHong Zhang     current_space_lvl->local_used      += nzk;
112178910aadSHong Zhang     current_space_lvl->local_remaining -= nzk;
112278910aadSHong Zhang 
112378910aadSHong Zhang     ui[k+1] = ui[k] + nzk;
112478910aadSHong Zhang   }
112578910aadSHong Zhang 
11269566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(perm,&rip));
11279566063dSJacob Faibussowitsch   PetscCall(PetscFree4(uj_ptr,uj_lvl_ptr,il,jl));
11289566063dSJacob Faibussowitsch   PetscCall(PetscFree(cols_lvl));
112978910aadSHong Zhang 
11309263d837SHong Zhang   /* copy free_space into uj and free free_space; set uj in new datastructure; */
11319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ui[am]+1,&uj));
11329566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceContiguous(&free_space,uj));
11339566063dSJacob Faibussowitsch   PetscCall(PetscIncompleteLLDestroy(lnk,lnkbt));
11349566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceDestroy(free_space_lvl));
113578910aadSHong Zhang 
113678910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1137719d5645SBarry Smith   B    = fact;
11389566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,NULL));
113978910aadSHong Zhang 
114078910aadSHong Zhang   b                = (Mat_SeqSBAIJ*)B->data;
114178910aadSHong Zhang   b->singlemalloc  = PETSC_FALSE;
1142e6b907acSBarry Smith   b->free_a        = PETSC_TRUE;
1143e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
114426fbe8dcSKarl Rupp 
11459566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ui[am]+1,&b->a));
114626fbe8dcSKarl Rupp 
114778910aadSHong Zhang   b->j             = uj;
114878910aadSHong Zhang   b->i             = ui;
1149f4259b30SLisandro Dalcin   b->diag          = NULL;
1150f4259b30SLisandro Dalcin   b->ilen          = NULL;
1151f4259b30SLisandro Dalcin   b->imax          = NULL;
115278910aadSHong Zhang   b->row           = perm;
115378910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
115426fbe8dcSKarl Rupp 
11559566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)perm));
115626fbe8dcSKarl Rupp 
115778910aadSHong Zhang   b->icol = perm;
115826fbe8dcSKarl Rupp 
11599566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)perm));
11609566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(am+1,&b->solve_work));
11619566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar))));
116226fbe8dcSKarl Rupp 
116378910aadSHong Zhang   b->maxnz = b->nz = ui[am];
116478910aadSHong Zhang 
116578910aadSHong Zhang   B->info.factor_mallocs   = reallocs;
116678910aadSHong Zhang   B->info.fill_ratio_given = fill;
116775567043SBarry Smith   if (ai[am] != 0.) {
116895b5ac22SHong Zhang     /* nonzeros in lower triangular part of A (includign diagonals)= (ai[am]+am)/2 */
116995b5ac22SHong Zhang     B->info.fill_ratio_needed = ((PetscReal)2*ui[am])/(ai[am]+am);
117078910aadSHong Zhang   } else {
117178910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
117278910aadSHong Zhang   }
11739263d837SHong Zhang #if defined(PETSC_USE_INFO)
11749263d837SHong Zhang   if (ai[am] != 0) {
117595b5ac22SHong Zhang     PetscReal af = B->info.fill_ratio_needed;
11769566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A,"Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af));
11779566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A,"Run with -pc_factor_fill %g or use \n",(double)af));
11789566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af));
11799263d837SHong Zhang   } else {
11809566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A,"Empty matrix\n"));
11819263d837SHong Zhang   }
11829263d837SHong Zhang #endif
118378910aadSHong Zhang   if (perm_identity) {
11840a3c351aSHong Zhang     B->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
11850a3c351aSHong Zhang     B->ops->solvetranspose        = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
118678910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
118778910aadSHong Zhang   } else {
1188719d5645SBarry Smith     (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
118978910aadSHong Zhang   }
1190c05c3958SHong Zhang   PetscFunctionReturn(0);
1191c05c3958SHong Zhang }
1192c05c3958SHong Zhang 
11930481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
1194c05c3958SHong Zhang {
119578910aadSHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
119678910aadSHong Zhang   Mat_SeqSBAIJ       *b;
119778910aadSHong Zhang   Mat                B;
11989186b105SHong Zhang   PetscBool          perm_identity,missing;
119978910aadSHong Zhang   PetscReal          fill = info->fill;
12005d0c19d7SBarry Smith   const PetscInt     *rip;
12015d0c19d7SBarry Smith   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,*ai=a->i,*aj=a->j,reallocs=0,prow;
120278910aadSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
120378910aadSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
12040298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
120578910aadSHong Zhang   PetscBT            lnkbt;
120678910aadSHong Zhang 
1207c05c3958SHong Zhang   PetscFunctionBegin;
12086ad2eaddSHong Zhang   if (bs > 1) { /* convert to seqsbaij */
12096ad2eaddSHong Zhang     if (!a->sbaijMat) {
12109566063dSJacob Faibussowitsch       PetscCall(MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat));
12116ad2eaddSHong Zhang     }
1212719d5645SBarry Smith     (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */
121326fbe8dcSKarl Rupp 
12149566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorSymbolic(fact,a->sbaijMat,perm,info));
12156ad2eaddSHong Zhang     PetscFunctionReturn(0);
12166ad2eaddSHong Zhang   }
12176ad2eaddSHong Zhang 
12189566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal(A,&missing,&i));
121928b400f6SJacob Faibussowitsch   PetscCheck(!missing,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %" PetscInt_FMT,i);
12209186b105SHong Zhang 
122178910aadSHong Zhang   /* check whether perm is the identity mapping */
12229566063dSJacob Faibussowitsch   PetscCall(ISIdentity(perm,&perm_identity));
122328b400f6SJacob Faibussowitsch   PetscCheck(perm_identity,PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported");
12249566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(perm,&rip));
122578910aadSHong Zhang 
122678910aadSHong Zhang   /* initialization */
12279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs+1,&ui));
1228e60cf9a0SBarry Smith   ui[0] = 0;
122978910aadSHong Zhang 
123078910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
123178910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
12329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(mbs,&ui_ptr,mbs,&il,mbs,&jl,mbs,&cols));
123378910aadSHong Zhang   for (i=0; i<mbs; i++) {
1234e60cf9a0SBarry Smith     jl[i] = mbs; il[i] = 0;
123578910aadSHong Zhang   }
123678910aadSHong Zhang 
123778910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
123878910aadSHong Zhang   nlnk = mbs + 1;
12399566063dSJacob Faibussowitsch   PetscCall(PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt));
124078910aadSHong Zhang 
12416a69fef8SHong Zhang   /* initial FreeSpace size is fill* (ai[mbs]+mbs)/2 */
12429566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[mbs]/2,mbs/2)),&free_space));
124326fbe8dcSKarl Rupp 
124478910aadSHong Zhang   current_space = free_space;
124578910aadSHong Zhang 
124678910aadSHong Zhang   for (k=0; k<mbs; k++) {  /* for each active row k */
124778910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
124878910aadSHong Zhang     nzk         = 0;
124978910aadSHong Zhang     ncols       = ai[rip[k]+1] - ai[rip[k]];
1250e60cf9a0SBarry Smith     ncols_upper = 0;
125178910aadSHong Zhang     for (j=0; j<ncols; j++) {
125278910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
125378910aadSHong Zhang       if (i >= k) { /* only take upper triangular entry */
125478910aadSHong Zhang         cols[ncols_upper] = i;
125578910aadSHong Zhang         ncols_upper++;
125678910aadSHong Zhang       }
125778910aadSHong Zhang     }
12589566063dSJacob Faibussowitsch     PetscCall(PetscLLAdd(ncols_upper,cols,mbs,&nlnk,lnk,lnkbt));
125978910aadSHong Zhang     nzk += nlnk;
126078910aadSHong Zhang 
126178910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
126278910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
126378910aadSHong Zhang 
126478910aadSHong Zhang     while (prow < k) {
126578910aadSHong Zhang       nextprow = jl[prow];
126678910aadSHong Zhang       /* merge prow into k-th row */
126778910aadSHong Zhang       jmin   = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
126878910aadSHong Zhang       jmax   = ui[prow+1];
126978910aadSHong Zhang       ncols  = jmax-jmin;
127078910aadSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
12719566063dSJacob Faibussowitsch       PetscCall(PetscLLAddSorted(ncols,uj_ptr,mbs,&nlnk,lnk,lnkbt));
127278910aadSHong Zhang       nzk   += nlnk;
127378910aadSHong Zhang 
127478910aadSHong Zhang       /* update il and jl for prow */
127578910aadSHong Zhang       if (jmin < jmax) {
127678910aadSHong Zhang         il[prow] = jmin;
127726fbe8dcSKarl Rupp         j        = *uj_ptr;
127826fbe8dcSKarl Rupp         jl[prow] = jl[j];
127926fbe8dcSKarl Rupp         jl[j]    = prow;
128078910aadSHong Zhang       }
128178910aadSHong Zhang       prow = nextprow;
128278910aadSHong Zhang     }
128378910aadSHong Zhang 
128478910aadSHong Zhang     /* if free space is not available, make more free space */
128578910aadSHong Zhang     if (current_space->local_remaining<nzk) {
128678910aadSHong Zhang       i    = mbs - k + 1; /* num of unfactored rows */
1287f91af8c7SBarry Smith       i    = PetscMin(PetscIntMultTruncate(i,nzk), PetscIntMultTruncate(i,i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
12889566063dSJacob Faibussowitsch       PetscCall(PetscFreeSpaceGet(i,&current_space));
128978910aadSHong Zhang       reallocs++;
129078910aadSHong Zhang     }
129178910aadSHong Zhang 
129278910aadSHong Zhang     /* copy data into free space, then initialize lnk */
12939566063dSJacob Faibussowitsch     PetscCall(PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt));
129478910aadSHong Zhang 
129578910aadSHong Zhang     /* add the k-th row into il and jl */
129678910aadSHong Zhang     if (nzk-1 > 0) {
129778910aadSHong Zhang       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
129878910aadSHong Zhang       jl[k] = jl[i]; jl[i] = k;
129978910aadSHong Zhang       il[k] = ui[k] + 1;
130078910aadSHong Zhang     }
130178910aadSHong Zhang     ui_ptr[k]                       = current_space->array;
130278910aadSHong Zhang     current_space->array           += nzk;
130378910aadSHong Zhang     current_space->local_used      += nzk;
130478910aadSHong Zhang     current_space->local_remaining -= nzk;
130578910aadSHong Zhang 
130678910aadSHong Zhang     ui[k+1] = ui[k] + nzk;
130778910aadSHong Zhang   }
130878910aadSHong Zhang 
13099566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(perm,&rip));
13109566063dSJacob Faibussowitsch   PetscCall(PetscFree4(ui_ptr,il,jl,cols));
131178910aadSHong Zhang 
13129263d837SHong Zhang   /* copy free_space into uj and free free_space; set uj in new datastructure; */
13139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ui[mbs]+1,&uj));
13149566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceContiguous(&free_space,uj));
13159566063dSJacob Faibussowitsch   PetscCall(PetscLLDestroy(lnk,lnkbt));
131678910aadSHong Zhang 
131778910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1318719d5645SBarry Smith   B    = fact;
13199566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,NULL));
132078910aadSHong Zhang 
132178910aadSHong Zhang   b               = (Mat_SeqSBAIJ*)B->data;
132278910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
1323e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1324e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
132526fbe8dcSKarl Rupp 
13269566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ui[mbs]+1,&b->a));
132726fbe8dcSKarl Rupp 
132878910aadSHong Zhang   b->j             = uj;
132978910aadSHong Zhang   b->i             = ui;
1330f4259b30SLisandro Dalcin   b->diag          = NULL;
1331f4259b30SLisandro Dalcin   b->ilen          = NULL;
1332f4259b30SLisandro Dalcin   b->imax          = NULL;
133378910aadSHong Zhang   b->row           = perm;
133478910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
133526fbe8dcSKarl Rupp 
13369566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)perm));
133778910aadSHong Zhang   b->icol  = perm;
13389566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)perm));
13399566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs+1,&b->solve_work));
13409566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar))));
134178910aadSHong Zhang   b->maxnz = b->nz = ui[mbs];
134278910aadSHong Zhang 
134378910aadSHong Zhang   B->info.factor_mallocs   = reallocs;
134478910aadSHong Zhang   B->info.fill_ratio_given = fill;
134575567043SBarry Smith   if (ai[mbs] != 0.) {
134695b5ac22SHong Zhang     /* nonzeros in lower triangular part of A = (ai[mbs]+mbs)/2 */
134795b5ac22SHong Zhang     B->info.fill_ratio_needed = ((PetscReal)2*ui[mbs])/(ai[mbs]+mbs);
134878910aadSHong Zhang   } else {
134978910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
135078910aadSHong Zhang   }
13519263d837SHong Zhang #if defined(PETSC_USE_INFO)
13529263d837SHong Zhang   if (ai[mbs] != 0.) {
13539263d837SHong Zhang     PetscReal af = B->info.fill_ratio_needed;
13549566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A,"Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af));
13559566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A,"Run with -pc_factor_fill %g or use \n",(double)af));
13569566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af));
13579263d837SHong Zhang   } else {
13589566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A,"Empty matrix\n"));
13599263d837SHong Zhang   }
13609263d837SHong Zhang #endif
136178910aadSHong Zhang   if (perm_identity) {
13626ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
136378910aadSHong Zhang   } else {
13646ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
136578910aadSHong Zhang   }
1366c05c3958SHong Zhang   PetscFunctionReturn(0);
1367c05c3958SHong Zhang }
1368c8342467SHong Zhang 
13694dd39f65SShri Abhyankar PetscErrorCode MatSolve_SeqBAIJ_N_NaturalOrdering(Mat A,Vec bb,Vec xx)
13701a83e813SShri Abhyankar {
13711a83e813SShri Abhyankar   Mat_SeqBAIJ       *a=(Mat_SeqBAIJ*)A->data;
13721a83e813SShri Abhyankar   const PetscInt    *ai=a->i,*aj=a->j,*adiag=a->diag,*vi;
13731a83e813SShri Abhyankar   PetscInt          i,k,n=a->mbs;
13741a83e813SShri Abhyankar   PetscInt          nz,bs=A->rmap->bs,bs2=a->bs2;
1375d9ca1df4SBarry Smith   const MatScalar   *aa=a->a,*v;
1376d9ca1df4SBarry Smith   PetscScalar       *x,*s,*t,*ls;
1377d9ca1df4SBarry Smith   const PetscScalar *b;
13781a83e813SShri Abhyankar 
13791a83e813SShri Abhyankar   PetscFunctionBegin;
13809566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb,&b));
13819566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx,&x));
13821a83e813SShri Abhyankar   t    = a->solve_work;
13831a83e813SShri Abhyankar 
13841a83e813SShri Abhyankar   /* forward solve the lower triangular */
13859566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t,b,bs)); /* copy 1st block of b to t */
13861a83e813SShri Abhyankar 
13871a83e813SShri Abhyankar   for (i=1; i<n; i++) {
13881a83e813SShri Abhyankar     v    = aa + bs2*ai[i];
13891a83e813SShri Abhyankar     vi   = aj + ai[i];
13901a83e813SShri Abhyankar     nz   = ai[i+1] - ai[i];
13911a83e813SShri Abhyankar     s    = t + bs*i;
13929566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(s,b+bs*i,bs)); /* copy i_th block of b to t */
13931a83e813SShri Abhyankar     for (k=0;k<nz;k++) {
139496b95a6bSBarry Smith       PetscKernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[k]);
13951a83e813SShri Abhyankar       v += bs2;
13961a83e813SShri Abhyankar     }
13971a83e813SShri Abhyankar   }
13981a83e813SShri Abhyankar 
13991a83e813SShri Abhyankar   /* backward solve the upper triangular */
14001a83e813SShri Abhyankar   ls = a->solve_work + A->cmap->n;
14011a83e813SShri Abhyankar   for (i=n-1; i>=0; i--) {
14021a83e813SShri Abhyankar     v    = aa + bs2*(adiag[i+1]+1);
14031a83e813SShri Abhyankar     vi   = aj + adiag[i+1]+1;
14041a83e813SShri Abhyankar     nz   = adiag[i] - adiag[i+1]-1;
14059566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(ls,t+i*bs,bs));
14061a83e813SShri Abhyankar     for (k=0; k<nz; k++) {
140796b95a6bSBarry Smith       PetscKernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[k]);
14081a83e813SShri Abhyankar       v += bs2;
14091a83e813SShri Abhyankar     }
141096b95a6bSBarry Smith     PetscKernel_w_gets_A_times_v(bs,ls,aa+bs2*adiag[i],t+i*bs); /* *inv(diagonal[i]) */
14119566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(x+i*bs,t+i*bs,bs));
14121a83e813SShri Abhyankar   }
14131a83e813SShri Abhyankar 
14149566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb,&b));
14159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx,&x));
14169566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n));
14171a83e813SShri Abhyankar   PetscFunctionReturn(0);
14181a83e813SShri Abhyankar }
14191a83e813SShri Abhyankar 
14204dd39f65SShri Abhyankar PetscErrorCode MatSolve_SeqBAIJ_N(Mat A,Vec bb,Vec xx)
142135aa4fcfSShri Abhyankar {
142235aa4fcfSShri Abhyankar   Mat_SeqBAIJ        *a   =(Mat_SeqBAIJ*)A->data;
142335aa4fcfSShri Abhyankar   IS                 iscol=a->col,isrow=a->row;
142435aa4fcfSShri Abhyankar   const PetscInt     *r,*c,*rout,*cout,*ai=a->i,*aj=a->j,*adiag=a->diag,*vi;
142535aa4fcfSShri Abhyankar   PetscInt           i,m,n=a->mbs;
142635aa4fcfSShri Abhyankar   PetscInt           nz,bs=A->rmap->bs,bs2=a->bs2;
1427d9ca1df4SBarry Smith   const MatScalar    *aa=a->a,*v;
1428d9ca1df4SBarry Smith   PetscScalar        *x,*s,*t,*ls;
1429d9ca1df4SBarry Smith   const PetscScalar  *b;
143035aa4fcfSShri Abhyankar 
143135aa4fcfSShri Abhyankar   PetscFunctionBegin;
14329566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb,&b));
14339566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx,&x));
143435aa4fcfSShri Abhyankar   t    = a->solve_work;
143535aa4fcfSShri Abhyankar 
14369566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow,&rout)); r = rout;
14379566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(iscol,&cout)); c = cout;
143835aa4fcfSShri Abhyankar 
143935aa4fcfSShri Abhyankar   /* forward solve the lower triangular */
14409566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t,b+bs*r[0],bs));
144135aa4fcfSShri Abhyankar   for (i=1; i<n; i++) {
144235aa4fcfSShri Abhyankar     v    = aa + bs2*ai[i];
144335aa4fcfSShri Abhyankar     vi   = aj + ai[i];
144435aa4fcfSShri Abhyankar     nz   = ai[i+1] - ai[i];
144535aa4fcfSShri Abhyankar     s    = t + bs*i;
14469566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(s,b+bs*r[i],bs));
144735aa4fcfSShri Abhyankar     for (m=0; m<nz; m++) {
144896b95a6bSBarry Smith       PetscKernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[m]);
144935aa4fcfSShri Abhyankar       v += bs2;
145035aa4fcfSShri Abhyankar     }
145135aa4fcfSShri Abhyankar   }
145235aa4fcfSShri Abhyankar 
145335aa4fcfSShri Abhyankar   /* backward solve the upper triangular */
145435aa4fcfSShri Abhyankar   ls = a->solve_work + A->cmap->n;
145535aa4fcfSShri Abhyankar   for (i=n-1; i>=0; i--) {
145635aa4fcfSShri Abhyankar     v    = aa + bs2*(adiag[i+1]+1);
145735aa4fcfSShri Abhyankar     vi   = aj + adiag[i+1]+1;
145835aa4fcfSShri Abhyankar     nz   = adiag[i] - adiag[i+1] - 1;
14599566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(ls,t+i*bs,bs));
146035aa4fcfSShri Abhyankar     for (m=0; m<nz; m++) {
146196b95a6bSBarry Smith       PetscKernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[m]);
146235aa4fcfSShri Abhyankar       v += bs2;
146335aa4fcfSShri Abhyankar     }
146496b95a6bSBarry Smith     PetscKernel_w_gets_A_times_v(bs,ls,v,t+i*bs); /* *inv(diagonal[i]) */
14659566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(x + bs*c[i],t+i*bs,bs));
146635aa4fcfSShri Abhyankar   }
14679566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow,&rout));
14689566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscol,&cout));
14699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb,&b));
14709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx,&x));
14719566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n));
147235aa4fcfSShri Abhyankar   PetscFunctionReturn(0);
147335aa4fcfSShri Abhyankar }
147435aa4fcfSShri Abhyankar 
1475a25532f0SBarry Smith /*
1476a25532f0SBarry Smith     For each block in an block array saves the largest absolute value in the block into another array
1477a25532f0SBarry Smith */
1478a25532f0SBarry Smith static PetscErrorCode MatBlockAbs_private(PetscInt nbs,PetscInt bs2,PetscScalar *blockarray,PetscReal *absarray)
14792efa7f71SHong Zhang {
14802efa7f71SHong Zhang   PetscInt       i,j;
14815fd66863SKarl Rupp 
14822efa7f71SHong Zhang   PetscFunctionBegin;
14839566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(absarray,nbs+1));
14842efa7f71SHong Zhang   for (i=0; i<nbs; i++) {
14852efa7f71SHong Zhang     for (j=0; j<bs2; j++) {
14862efa7f71SHong Zhang       if (absarray[i] < PetscAbsScalar(blockarray[i*nbs+j])) absarray[i] = PetscAbsScalar(blockarray[i*nbs+j]);
14872efa7f71SHong Zhang     }
14882efa7f71SHong Zhang   }
14892efa7f71SHong Zhang   PetscFunctionReturn(0);
14902efa7f71SHong Zhang }
14912efa7f71SHong Zhang 
1492fe97e370SBarry Smith /*
1493fe97e370SBarry Smith      This needs to be renamed and called by the regular MatILUFactor_SeqBAIJ when drop tolerance is used
1494fe97e370SBarry Smith */
1495c8342467SHong Zhang PetscErrorCode MatILUDTFactor_SeqBAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact)
1496c8342467SHong Zhang {
14972efa7f71SHong Zhang   Mat            B = *fact;
14982efa7f71SHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b;
14992efa7f71SHong Zhang   IS             isicol;
15002efa7f71SHong Zhang   const PetscInt *r,*ic;
15012efa7f71SHong Zhang   PetscInt       i,mbs=a->mbs,bs=A->rmap->bs,bs2=a->bs2,*ai=a->i,*aj=a->j,*ajtmp,*adiag;
15022efa7f71SHong Zhang   PetscInt       *bi,*bj,*bdiag;
15032efa7f71SHong Zhang 
15042efa7f71SHong Zhang   PetscInt  row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au;
15052efa7f71SHong Zhang   PetscInt  nlnk,*lnk;
15062efa7f71SHong Zhang   PetscBT   lnkbt;
1507ace3abfcSBarry Smith   PetscBool row_identity,icol_identity;
15082efa7f71SHong Zhang   MatScalar *aatmp,*pv,*batmp,*ba,*rtmp,*pc,*multiplier,*vtmp;
15092efa7f71SHong Zhang   PetscInt  j,nz,*pj,*bjtmp,k,ncut,*jtmp;
15102efa7f71SHong Zhang 
1511182b8fbaSHong Zhang   PetscReal dt=info->dt;          /* shift=info->shiftamount; */
15122efa7f71SHong Zhang   PetscInt  nnz_max;
1513ace3abfcSBarry Smith   PetscBool missing;
1514021822bcSHong Zhang   PetscReal *vtmp_abs;
151597ef94ebSSatish Balay   MatScalar *v_work;
151697ef94ebSSatish Balay   PetscInt  *v_pivots;
15175f8bbccaSHong Zhang   PetscBool allowzeropivot,zeropivotdetected=PETSC_FALSE;
15182efa7f71SHong Zhang 
1519c8342467SHong Zhang   PetscFunctionBegin;
15202efa7f71SHong Zhang   /* ------- symbolic factorization, can be reused ---------*/
15219566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal(A,&missing,&i));
152228b400f6SJacob Faibussowitsch   PetscCheck(!missing,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %" PetscInt_FMT,i);
15232efa7f71SHong Zhang   adiag=a->diag;
15242efa7f71SHong Zhang 
15259566063dSJacob Faibussowitsch   PetscCall(ISInvertPermutation(iscol,PETSC_DECIDE,&isicol));
15262efa7f71SHong Zhang 
15272efa7f71SHong Zhang   /* bdiag is location of diagonal in factor */
15289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs+1,&bdiag));
15292efa7f71SHong Zhang 
15302efa7f71SHong Zhang   /* allocate row pointers bi */
15319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2*mbs+2,&bi));
15322efa7f71SHong Zhang 
15332efa7f71SHong Zhang   /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */
15342efa7f71SHong Zhang   dtcount = (PetscInt)info->dtcount;
15356bce7ff8SHong Zhang   if (dtcount > mbs-1) dtcount = mbs-1;
15366bce7ff8SHong Zhang   nnz_max = ai[mbs]+2*mbs*dtcount +2;
15376da515a1SHong Zhang   /* printf("MatILUDTFactor_SeqBAIJ, bs %d, ai[mbs] %d, nnz_max  %d, dtcount %d\n",bs,ai[mbs],nnz_max,dtcount); */
15389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nnz_max,&bj));
15396bce7ff8SHong Zhang   nnz_max = nnz_max*bs2;
15409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nnz_max,&ba));
15412efa7f71SHong Zhang 
15422efa7f71SHong Zhang   /* put together the new matrix */
15439566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,NULL));
15449566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)B,(PetscObject)isicol));
154526fbe8dcSKarl Rupp 
15462efa7f71SHong Zhang   b               = (Mat_SeqBAIJ*)(B)->data;
15472efa7f71SHong Zhang   b->free_a       = PETSC_TRUE;
15482efa7f71SHong Zhang   b->free_ij      = PETSC_TRUE;
15492efa7f71SHong Zhang   b->singlemalloc = PETSC_FALSE;
155026fbe8dcSKarl Rupp 
15512efa7f71SHong Zhang   b->a    = ba;
15522efa7f71SHong Zhang   b->j    = bj;
15532efa7f71SHong Zhang   b->i    = bi;
15542efa7f71SHong Zhang   b->diag = bdiag;
1555f4259b30SLisandro Dalcin   b->ilen = NULL;
1556f4259b30SLisandro Dalcin   b->imax = NULL;
15572efa7f71SHong Zhang   b->row  = isrow;
15582efa7f71SHong Zhang   b->col  = iscol;
155926fbe8dcSKarl Rupp 
15609566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)isrow));
15619566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)iscol));
156226fbe8dcSKarl Rupp 
15632efa7f71SHong Zhang   b->icol  = isicol;
15649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs*(mbs+1),&b->solve_work));
15659566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar))));
15662efa7f71SHong Zhang   b->maxnz = nnz_max/bs2;
15672efa7f71SHong Zhang 
1568d5f3da31SBarry Smith   (B)->factortype            = MAT_FACTOR_ILUDT;
15692efa7f71SHong Zhang   (B)->info.factor_mallocs   = 0;
15702efa7f71SHong Zhang   (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)(ai[mbs]*bs2));
15712efa7f71SHong Zhang   /* ------- end of symbolic factorization ---------*/
15729566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow,&r));
15739566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol,&ic));
15742efa7f71SHong Zhang 
15752efa7f71SHong Zhang   /* linked list for storing column indices of the active row */
15762efa7f71SHong Zhang   nlnk = mbs + 1;
15779566063dSJacob Faibussowitsch   PetscCall(PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt));
15782efa7f71SHong Zhang 
15792efa7f71SHong Zhang   /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */
15809566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(mbs,&im,mbs,&jtmp));
15812efa7f71SHong Zhang   /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */
15829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(mbs*bs2,&rtmp,mbs*bs2,&vtmp));
15839566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs+1,&vtmp_abs));
15849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(bs,&v_work,bs2,&multiplier,bs,&v_pivots));
15852efa7f71SHong Zhang 
15865f8bbccaSHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
1587e60cf9a0SBarry Smith   bi[0]       = 0;
15882efa7f71SHong Zhang   bdiag[0]    = (nnz_max/bs2)-1; /* location of diagonal in factor B */
15896bce7ff8SHong Zhang   bi[2*mbs+1] = bdiag[0]+1; /* endof bj and ba array */
15902efa7f71SHong Zhang   for (i=0; i<mbs; i++) {
15912efa7f71SHong Zhang     /* copy initial fill into linked list */
15922efa7f71SHong Zhang     nzi = ai[r[i]+1] - ai[r[i]];
159328b400f6SJacob 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);
15942efa7f71SHong Zhang     nzi_al = adiag[r[i]] - ai[r[i]];
15952efa7f71SHong Zhang     nzi_au = ai[r[i]+1] - adiag[r[i]] -1;
15962efa7f71SHong Zhang 
15972efa7f71SHong Zhang     /* load in initial unfactored row */
15982efa7f71SHong Zhang     ajtmp = aj + ai[r[i]];
15999566063dSJacob Faibussowitsch     PetscCall(PetscLLAddPerm(nzi,ajtmp,ic,mbs,&nlnk,lnk,lnkbt));
16009566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(rtmp,mbs*bs2));
16012efa7f71SHong Zhang     aatmp = a->a + bs2*ai[r[i]];
16029566063dSJacob Faibussowitsch     for (j=0; j<nzi; j++) PetscCall(PetscArraycpy(rtmp+bs2*ic[ajtmp[j]],aatmp+bs2*j,bs2));
16032efa7f71SHong Zhang 
16042efa7f71SHong Zhang     /* add pivot rows into linked list */
16052efa7f71SHong Zhang     row = lnk[mbs];
16062efa7f71SHong Zhang     while (row < i) {
16072efa7f71SHong Zhang       nzi_bl = bi[row+1] - bi[row] + 1;
16082efa7f71SHong Zhang       bjtmp  = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */
16099566063dSJacob Faibussowitsch       PetscCall(PetscLLAddSortedLU(bjtmp,row,&nlnk,lnk,lnkbt,i,nzi_bl,im));
16102efa7f71SHong Zhang       nzi   += nlnk;
16112efa7f71SHong Zhang       row    = lnk[row];
16122efa7f71SHong Zhang     }
16132efa7f71SHong Zhang 
16142efa7f71SHong Zhang     /* copy data from lnk into jtmp, then initialize lnk */
16159566063dSJacob Faibussowitsch     PetscCall(PetscLLClean(mbs,mbs,nzi,lnk,jtmp,lnkbt));
16162efa7f71SHong Zhang 
16172efa7f71SHong Zhang     /* numerical factorization */
16182efa7f71SHong Zhang     bjtmp = jtmp;
16192efa7f71SHong Zhang     row   = *bjtmp++; /* 1st pivot row */
16202efa7f71SHong Zhang 
16212efa7f71SHong Zhang     while  (row < i) {
16222efa7f71SHong Zhang       pc = rtmp + bs2*row;
16232efa7f71SHong Zhang       pv = ba + bs2*bdiag[row]; /* inv(diag) of the pivot row */
162496b95a6bSBarry Smith       PetscKernel_A_gets_A_times_B(bs,pc,pv,multiplier); /* pc= multiplier = pc*inv(diag[row]) */
16259566063dSJacob Faibussowitsch       PetscCall(MatBlockAbs_private(1,bs2,pc,vtmp_abs));
16262efa7f71SHong Zhang       if (vtmp_abs[0] > dt) { /* apply tolerance dropping rule */
16272efa7f71SHong Zhang         pj = bj + bdiag[row+1] + 1;         /* point to 1st entry of U(row,:) */
16282efa7f71SHong Zhang         pv = ba + bs2*(bdiag[row+1] + 1);
16292efa7f71SHong Zhang         nz = bdiag[row] - bdiag[row+1] - 1;         /* num of entries in U(row,:), excluding diagonal */
16302efa7f71SHong Zhang         for (j=0; j<nz; j++) {
163196b95a6bSBarry Smith           PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j);
16322efa7f71SHong Zhang         }
16339566063dSJacob Faibussowitsch         /* PetscCall(PetscLogFlops(bslog*(nz+1.0)-bs)); */
16342efa7f71SHong Zhang       }
16352efa7f71SHong Zhang       row = *bjtmp++;
16362efa7f71SHong Zhang     }
16372efa7f71SHong Zhang 
16382efa7f71SHong Zhang     /* copy sparse rtmp into contiguous vtmp; separate L and U part */
16392efa7f71SHong Zhang     nzi_bl = 0; j = 0;
16402efa7f71SHong Zhang     while (jtmp[j] < i) { /* L-part. Note: jtmp is sorted */
16419566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2));
16422efa7f71SHong Zhang       nzi_bl++; j++;
16432efa7f71SHong Zhang     }
16442efa7f71SHong Zhang     nzi_bu = nzi - nzi_bl -1;
16452efa7f71SHong Zhang 
16462efa7f71SHong Zhang     while (j < nzi) { /* U-part */
16479566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2));
16482efa7f71SHong Zhang       j++;
16492efa7f71SHong Zhang     }
16502efa7f71SHong Zhang 
16519566063dSJacob Faibussowitsch     PetscCall(MatBlockAbs_private(nzi,bs2,vtmp,vtmp_abs));
16525f8bbccaSHong Zhang 
16532efa7f71SHong Zhang     bjtmp = bj + bi[i];
16542efa7f71SHong Zhang     batmp = ba + bs2*bi[i];
16552efa7f71SHong Zhang     /* apply level dropping rule to L part */
16562efa7f71SHong Zhang     ncut = nzi_al + dtcount;
16572efa7f71SHong Zhang     if (ncut < nzi_bl) {
16589566063dSJacob Faibussowitsch       PetscCall(PetscSortSplitReal(ncut,nzi_bl,vtmp_abs,jtmp));
16599566063dSJacob Faibussowitsch       PetscCall(PetscSortIntWithScalarArray(ncut,jtmp,vtmp));
16602efa7f71SHong Zhang     } else {
16612efa7f71SHong Zhang       ncut = nzi_bl;
16622efa7f71SHong Zhang     }
16632efa7f71SHong Zhang     for (j=0; j<ncut; j++) {
16642efa7f71SHong Zhang       bjtmp[j] = jtmp[j];
16659566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(batmp+bs2*j,rtmp+bs2*bjtmp[j],bs2));
16662efa7f71SHong Zhang     }
16672efa7f71SHong Zhang     bi[i+1] = bi[i] + ncut;
16682efa7f71SHong Zhang     nzi     = ncut + 1;
16692efa7f71SHong Zhang 
16702efa7f71SHong Zhang     /* apply level dropping rule to U part */
16712efa7f71SHong Zhang     ncut = nzi_au + dtcount;
16722efa7f71SHong Zhang     if (ncut < nzi_bu) {
16739566063dSJacob Faibussowitsch       PetscCall(PetscSortSplitReal(ncut,nzi_bu,vtmp_abs+nzi_bl+1,jtmp+nzi_bl+1));
16749566063dSJacob Faibussowitsch       PetscCall(PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1));
16752efa7f71SHong Zhang     } else {
16762efa7f71SHong Zhang       ncut = nzi_bu;
16772efa7f71SHong Zhang     }
16782efa7f71SHong Zhang     nzi += ncut;
16792efa7f71SHong Zhang 
16802efa7f71SHong Zhang     /* mark bdiagonal */
16812efa7f71SHong Zhang     bdiag[i+1]    = bdiag[i] - (ncut + 1);
16826bce7ff8SHong Zhang     bi[2*mbs - i] = bi[2*mbs - i +1] - (ncut + 1);
16836bce7ff8SHong Zhang 
16842efa7f71SHong Zhang     bjtmp  = bj + bdiag[i];
16852efa7f71SHong Zhang     batmp  = ba + bs2*bdiag[i];
16869566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(batmp,rtmp+bs2*i,bs2));
16872efa7f71SHong Zhang     *bjtmp = i;
16885f8bbccaSHong Zhang 
16892efa7f71SHong Zhang     bjtmp = bj + bdiag[i+1]+1;
16902efa7f71SHong Zhang     batmp = ba + (bdiag[i+1]+1)*bs2;
16912efa7f71SHong Zhang 
16922efa7f71SHong Zhang     for (k=0; k<ncut; k++) {
16932efa7f71SHong Zhang       bjtmp[k] = jtmp[nzi_bl+1+k];
16949566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(batmp+bs2*k,rtmp+bs2*bjtmp[k],bs2));
16952efa7f71SHong Zhang     }
16962efa7f71SHong Zhang 
16972efa7f71SHong Zhang     im[i] = nzi; /* used by PetscLLAddSortedLU() */
16982efa7f71SHong Zhang 
1699a5b23f4aSJose E. Roman     /* invert diagonal block for simpler triangular solves - add shift??? */
17002efa7f71SHong Zhang     batmp = ba + bs2*bdiag[i];
17015f8bbccaSHong Zhang 
17029566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A(bs,batmp,v_pivots,v_work,allowzeropivot,&zeropivotdetected));
17037b6c816cSBarry Smith     if (zeropivotdetected) B->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
17042efa7f71SHong Zhang   } /* for (i=0; i<mbs; i++) */
17059566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v_work,multiplier,v_pivots));
17062efa7f71SHong Zhang 
17076da515a1SHong Zhang   /* printf("end of L %d, beginning of U %d\n",bi[mbs],bdiag[mbs]); */
170894bad497SJacob 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]);
17092efa7f71SHong Zhang 
17109566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow,&r));
17119566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol,&ic));
17122efa7f71SHong Zhang 
17139566063dSJacob Faibussowitsch   PetscCall(PetscLLDestroy(lnk,lnkbt));
17142efa7f71SHong Zhang 
17159566063dSJacob Faibussowitsch   PetscCall(PetscFree2(im,jtmp));
17169566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rtmp,vtmp));
17172efa7f71SHong Zhang 
17189566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(bs2*B->cmap->n));
17192efa7f71SHong Zhang   b->maxnz = b->nz = bi[mbs] + bdiag[0] - bdiag[mbs];
17202efa7f71SHong Zhang 
17219566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isrow,&row_identity));
17229566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isicol,&icol_identity));
17232efa7f71SHong Zhang   if (row_identity && icol_identity) {
17244dd39f65SShri Abhyankar     B->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering;
17252efa7f71SHong Zhang   } else {
17264dd39f65SShri Abhyankar     B->ops->solve = MatSolve_SeqBAIJ_N;
17272efa7f71SHong Zhang   }
17282efa7f71SHong Zhang 
1729f4259b30SLisandro Dalcin   B->ops->solveadd          = NULL;
1730f4259b30SLisandro Dalcin   B->ops->solvetranspose    = NULL;
1731f4259b30SLisandro Dalcin   B->ops->solvetransposeadd = NULL;
1732f4259b30SLisandro Dalcin   B->ops->matsolve          = NULL;
17332efa7f71SHong Zhang   B->assembled              = PETSC_TRUE;
17342efa7f71SHong Zhang   B->preallocated           = PETSC_TRUE;
1735c8342467SHong Zhang   PetscFunctionReturn(0);
1736c8342467SHong Zhang }
1737