1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 35d517e7eSBarry Smith /* 4ec3a8b7bSBarry Smith Factorization code for BAIJ format. 55d517e7eSBarry Smith */ 67c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h" 7c60f0209SBarry Smith #include "../src/mat/blockinvert.h" 8ec3a8b7bSBarry Smith 9b588c5a2SHong Zhang #undef __FUNCT__ 104dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2" 114dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat B,Mat A,const MatFactorInfo *info) 12b588c5a2SHong Zhang { 13b588c5a2SHong Zhang Mat C=B; 14b588c5a2SHong Zhang Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; 15b588c5a2SHong Zhang IS isrow = b->row,isicol = b->icol; 16b588c5a2SHong Zhang PetscErrorCode ierr; 17b588c5a2SHong Zhang const PetscInt *r,*ic,*ics; 18bbd65245SShri Abhyankar PetscInt i,j,k,nz,nzL,row,*pj; 19bbd65245SShri Abhyankar const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,bs2=a->bs2; 20bbd65245SShri Abhyankar const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag; 21bbd65245SShri Abhyankar MatScalar *rtmp,*pc,*mwork,*pv; 22bbd65245SShri Abhyankar MatScalar *aa=a->a,*v; 23bbd65245SShri Abhyankar PetscInt flg; 24182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 25b588c5a2SHong Zhang 26b588c5a2SHong Zhang PetscFunctionBegin; 27b588c5a2SHong Zhang ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 28b588c5a2SHong Zhang ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 29b588c5a2SHong Zhang 304c000e2eSHong Zhang /* generate work space needed by the factorization */ 31fca92195SBarry Smith ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr); 324c000e2eSHong Zhang ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 33b588c5a2SHong Zhang ics = ic; 34b588c5a2SHong Zhang 35b588c5a2SHong Zhang for (i=0; i<n; i++){ 36b588c5a2SHong Zhang /* zero rtmp */ 37b588c5a2SHong Zhang /* L part */ 38b588c5a2SHong Zhang nz = bi[i+1] - bi[i]; 39b588c5a2SHong Zhang bjtmp = bj + bi[i]; 40b588c5a2SHong Zhang for (j=0; j<nz; j++){ 41b588c5a2SHong Zhang ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 42b588c5a2SHong Zhang } 43b588c5a2SHong Zhang 44b588c5a2SHong Zhang /* U part */ 450c4413a7SShri Abhyankar nz = bdiag[i] - bdiag[i+1]; 460c4413a7SShri Abhyankar bjtmp = bj + bdiag[i+1]+1; 470c4413a7SShri Abhyankar for (j=0; j<nz; j++){ 480c4413a7SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 490c4413a7SShri Abhyankar } 500c4413a7SShri Abhyankar 510c4413a7SShri Abhyankar /* load in initial (unfactored row) */ 520c4413a7SShri Abhyankar nz = ai[r[i]+1] - ai[r[i]]; 530c4413a7SShri Abhyankar ajtmp = aj + ai[r[i]]; 540c4413a7SShri Abhyankar v = aa + bs2*ai[r[i]]; 550c4413a7SShri Abhyankar for (j=0; j<nz; j++) { 560c4413a7SShri Abhyankar ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 570c4413a7SShri Abhyankar } 580c4413a7SShri Abhyankar 590c4413a7SShri Abhyankar /* elimination */ 600c4413a7SShri Abhyankar bjtmp = bj + bi[i]; 610c4413a7SShri Abhyankar nzL = bi[i+1] - bi[i]; 620c4413a7SShri Abhyankar for(k=0;k < nzL;k++) { 630c4413a7SShri Abhyankar row = bjtmp[k]; 640c4413a7SShri Abhyankar pc = rtmp + bs2*row; 650c4413a7SShri Abhyankar for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }} 660c4413a7SShri Abhyankar if (flg) { 670c4413a7SShri Abhyankar pv = b->a + bs2*bdiag[row]; 680c4413a7SShri Abhyankar /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 690c4413a7SShri Abhyankar ierr = Kernel_A_gets_A_times_B_2(pc,pv,mwork);CHKERRQ(ierr); 700c4413a7SShri Abhyankar 710c4413a7SShri Abhyankar pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */ 720c4413a7SShri Abhyankar pv = b->a + bs2*(bdiag[row+1]+1); 730c4413a7SShri Abhyankar nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */ 740c4413a7SShri Abhyankar for (j=0; j<nz; j++) { 750c4413a7SShri Abhyankar /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 760c4413a7SShri Abhyankar /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 770c4413a7SShri Abhyankar v = rtmp + 4*pj[j]; 780c4413a7SShri Abhyankar ierr = Kernel_A_gets_A_minus_B_times_C_2(v,pc,pv);CHKERRQ(ierr); 790c4413a7SShri Abhyankar pv += 4; 800c4413a7SShri Abhyankar } 810c4413a7SShri Abhyankar ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 820c4413a7SShri Abhyankar } 830c4413a7SShri Abhyankar } 840c4413a7SShri Abhyankar 850c4413a7SShri Abhyankar /* finished row so stick it into b->a */ 860c4413a7SShri Abhyankar /* L part */ 870c4413a7SShri Abhyankar pv = b->a + bs2*bi[i] ; 880c4413a7SShri Abhyankar pj = b->j + bi[i] ; 890c4413a7SShri Abhyankar nz = bi[i+1] - bi[i]; 900c4413a7SShri Abhyankar for (j=0; j<nz; j++) { 910c4413a7SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 920c4413a7SShri Abhyankar } 930c4413a7SShri Abhyankar 940c4413a7SShri Abhyankar /* Mark diagonal and invert diagonal for simplier triangular solves */ 950c4413a7SShri Abhyankar pv = b->a + bs2*bdiag[i]; 960c4413a7SShri Abhyankar pj = b->j + bdiag[i]; 970c4413a7SShri Abhyankar ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 980c4413a7SShri Abhyankar /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 990c4413a7SShri Abhyankar ierr = Kernel_A_gets_inverse_A_2(pv,shift);CHKERRQ(ierr); 1000c4413a7SShri Abhyankar 1010c4413a7SShri Abhyankar /* U part */ 1020c4413a7SShri Abhyankar pv = b->a + bs2*(bdiag[i+1]+1); 1030c4413a7SShri Abhyankar pj = b->j + bdiag[i+1]+1; 1040c4413a7SShri Abhyankar nz = bdiag[i] - bdiag[i+1] - 1; 1050c4413a7SShri Abhyankar for (j=0; j<nz; j++){ 1060c4413a7SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1070c4413a7SShri Abhyankar } 1080c4413a7SShri Abhyankar } 1090c4413a7SShri Abhyankar 110fca92195SBarry Smith ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); 1110c4413a7SShri Abhyankar ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 1120c4413a7SShri Abhyankar ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 1134dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_2; 1144dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2; 1150c4413a7SShri Abhyankar 1160c4413a7SShri Abhyankar C->assembled = PETSC_TRUE; 117766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*2*2*2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 1180c4413a7SShri Abhyankar PetscFunctionReturn(0); 1190c4413a7SShri Abhyankar } 1200c4413a7SShri Abhyankar 1214c000e2eSHong Zhang #undef __FUNCT__ 1224dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering" 1234dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info) 1244c000e2eSHong Zhang { 1254c000e2eSHong Zhang Mat C=B; 1264c000e2eSHong Zhang Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; 1274c000e2eSHong Zhang PetscErrorCode ierr; 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; 1354c000e2eSHong Zhang 1364c000e2eSHong Zhang PetscFunctionBegin; 1374c000e2eSHong Zhang /* generate work space needed by the factorization */ 138fca92195SBarry Smith ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr); 1394c000e2eSHong Zhang ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 1404c000e2eSHong Zhang 1414c000e2eSHong Zhang for (i=0; i<n; i++){ 1424c000e2eSHong Zhang /* zero rtmp */ 1434c000e2eSHong Zhang /* L part */ 1444c000e2eSHong Zhang nz = bi[i+1] - bi[i]; 1454c000e2eSHong Zhang bjtmp = bj + bi[i]; 1464c000e2eSHong Zhang for (j=0; j<nz; j++){ 1474c000e2eSHong Zhang ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1484c000e2eSHong Zhang } 1494c000e2eSHong Zhang 1504c000e2eSHong Zhang /* U part */ 151b2b2dd24SShri Abhyankar nz = bdiag[i] - bdiag[i+1]; 152b2b2dd24SShri Abhyankar bjtmp = bj + bdiag[i+1]+1; 153b2b2dd24SShri Abhyankar for (j=0; j<nz; j++){ 154b2b2dd24SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 155b2b2dd24SShri Abhyankar } 156b2b2dd24SShri Abhyankar 157b2b2dd24SShri Abhyankar /* load in initial (unfactored row) */ 158b2b2dd24SShri Abhyankar nz = ai[i+1] - ai[i]; 159b2b2dd24SShri Abhyankar ajtmp = aj + ai[i]; 160b2b2dd24SShri Abhyankar v = aa + bs2*ai[i]; 161b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 162b2b2dd24SShri Abhyankar ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 163b2b2dd24SShri Abhyankar } 164b2b2dd24SShri Abhyankar 165b2b2dd24SShri Abhyankar /* elimination */ 166b2b2dd24SShri Abhyankar bjtmp = bj + bi[i]; 167b2b2dd24SShri Abhyankar nzL = bi[i+1] - bi[i]; 168b2b2dd24SShri Abhyankar for(k=0;k < nzL;k++) { 169b2b2dd24SShri Abhyankar row = bjtmp[k]; 170b2b2dd24SShri Abhyankar pc = rtmp + bs2*row; 171b2b2dd24SShri Abhyankar for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }} 172b2b2dd24SShri Abhyankar if (flg) { 173b2b2dd24SShri Abhyankar pv = b->a + bs2*bdiag[row]; 174b2b2dd24SShri Abhyankar /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 175b2b2dd24SShri Abhyankar ierr = Kernel_A_gets_A_times_B_2(pc,pv,mwork);CHKERRQ(ierr); 176b2b2dd24SShri Abhyankar 177b2b2dd24SShri Abhyankar pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */ 178b2b2dd24SShri Abhyankar pv = b->a + bs2*(bdiag[row+1]+1); 179b2b2dd24SShri Abhyankar nz = bdiag[row]-bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */ 180b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 181b2b2dd24SShri Abhyankar /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 182b2b2dd24SShri Abhyankar /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 183b2b2dd24SShri Abhyankar v = rtmp + 4*pj[j]; 184b2b2dd24SShri Abhyankar ierr = Kernel_A_gets_A_minus_B_times_C_2(v,pc,pv);CHKERRQ(ierr); 185b2b2dd24SShri Abhyankar pv += 4; 186b2b2dd24SShri Abhyankar } 187b2b2dd24SShri Abhyankar ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 188b2b2dd24SShri Abhyankar } 189b2b2dd24SShri Abhyankar } 190b2b2dd24SShri Abhyankar 191b2b2dd24SShri Abhyankar /* finished row so stick it into b->a */ 192b2b2dd24SShri Abhyankar /* L part */ 193b2b2dd24SShri Abhyankar pv = b->a + bs2*bi[i] ; 194b2b2dd24SShri Abhyankar pj = b->j + bi[i] ; 195b2b2dd24SShri Abhyankar nz = bi[i+1] - bi[i]; 196b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 197b2b2dd24SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 198b2b2dd24SShri Abhyankar } 199b2b2dd24SShri Abhyankar 200b2b2dd24SShri Abhyankar /* Mark diagonal and invert diagonal for simplier triangular solves */ 201b2b2dd24SShri Abhyankar pv = b->a + bs2*bdiag[i]; 202b2b2dd24SShri Abhyankar pj = b->j + bdiag[i]; 203b2b2dd24SShri Abhyankar ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 204b2b2dd24SShri Abhyankar /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 205b2b2dd24SShri Abhyankar ierr = Kernel_A_gets_inverse_A_2(pv,shift);CHKERRQ(ierr); 206b2b2dd24SShri Abhyankar 207b2b2dd24SShri Abhyankar /* U part */ 208b2b2dd24SShri Abhyankar /* 209b2b2dd24SShri Abhyankar pv = b->a + bs2*bi[2*n-i]; 210b2b2dd24SShri Abhyankar pj = b->j + bi[2*n-i]; 211b2b2dd24SShri Abhyankar nz = bi[2*n-i+1] - bi[2*n-i] - 1; 212b2b2dd24SShri Abhyankar */ 213b2b2dd24SShri Abhyankar pv = b->a + bs2*(bdiag[i+1]+1); 214b2b2dd24SShri Abhyankar pj = b->j + bdiag[i+1]+1; 215b2b2dd24SShri Abhyankar nz = bdiag[i] - bdiag[i+1] - 1; 216b2b2dd24SShri Abhyankar for (j=0; j<nz; j++){ 217b2b2dd24SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 218b2b2dd24SShri Abhyankar } 219b2b2dd24SShri Abhyankar } 220fca92195SBarry Smith ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); 221ae3d28f0SHong Zhang 2224dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering; 2234dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering; 224b2b2dd24SShri Abhyankar C->assembled = PETSC_TRUE; 225766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*2*2*2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 226b2b2dd24SShri Abhyankar PetscFunctionReturn(0); 227b2b2dd24SShri Abhyankar } 228b2b2dd24SShri Abhyankar 229b2b2dd24SShri Abhyankar #undef __FUNCT__ 23006e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_inplace" 23106e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_inplace(Mat B,Mat A,const MatFactorInfo *info) 2324eeb42bcSBarry Smith { 233719d5645SBarry Smith Mat C = B; 2344eeb42bcSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 2357cdf2dbbSSatish Balay IS isrow = b->row,isicol = b->icol; 2366849ba73SBarry Smith PetscErrorCode ierr; 2375d0c19d7SBarry Smith const PetscInt *r,*ic; 2385d0c19d7SBarry Smith PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 239690b6cddSBarry Smith PetscInt *ajtmpold,*ajtmp,nz,row; 240690b6cddSBarry Smith PetscInt *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj; 241329f5518SBarry Smith MatScalar *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4; 2422515f8d2SSatish Balay MatScalar p1,p2,p3,p4; 2433f1db9ecSBarry Smith MatScalar *ba = b->a,*aa = a->a; 244182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 2454eeb42bcSBarry Smith 2463a40ed3dSBarry Smith PetscFunctionBegin; 2474eeb42bcSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 2484eeb42bcSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 249b0a32e0cSBarry Smith ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 2504eeb42bcSBarry Smith 2514eeb42bcSBarry Smith for (i=0; i<n; i++) { 2524078e994SBarry Smith nz = bi[i+1] - bi[i]; 2534078e994SBarry Smith ajtmp = bj + bi[i]; 2544eeb42bcSBarry Smith for (j=0; j<nz; j++) { 2554eeb42bcSBarry Smith x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0; 2564eeb42bcSBarry Smith } 2574eeb42bcSBarry Smith /* load in initial (unfactored row) */ 2584eeb42bcSBarry Smith idx = r[i]; 2594078e994SBarry Smith nz = ai[idx+1] - ai[idx]; 2604078e994SBarry Smith ajtmpold = aj + ai[idx]; 2614078e994SBarry Smith v = aa + 4*ai[idx]; 2624eeb42bcSBarry Smith for (j=0; j<nz; j++) { 2634eeb42bcSBarry Smith x = rtmp+4*ic[ajtmpold[j]]; 2644eeb42bcSBarry Smith x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 2654eeb42bcSBarry Smith v += 4; 2664eeb42bcSBarry Smith } 2674eeb42bcSBarry Smith row = *ajtmp++; 2684eeb42bcSBarry Smith while (row < i) { 2694eeb42bcSBarry Smith pc = rtmp + 4*row; 2704eeb42bcSBarry Smith p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 27188685aaeSLois Curfman McInnes if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) { 2724078e994SBarry Smith pv = ba + 4*diag_offset[row]; 2734078e994SBarry Smith pj = bj + diag_offset[row] + 1; 2744eeb42bcSBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 2754eeb42bcSBarry Smith pc[0] = m1 = p1*x1 + p3*x2; 2764eeb42bcSBarry Smith pc[1] = m2 = p2*x1 + p4*x2; 2774eeb42bcSBarry Smith pc[2] = m3 = p1*x3 + p3*x4; 2784eeb42bcSBarry Smith pc[3] = m4 = p2*x3 + p4*x4; 2794078e994SBarry Smith nz = bi[row+1] - diag_offset[row] - 1; 2804eeb42bcSBarry Smith pv += 4; 2814eeb42bcSBarry Smith for (j=0; j<nz; j++) { 2824eeb42bcSBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 2834eeb42bcSBarry Smith x = rtmp + 4*pj[j]; 2844eeb42bcSBarry Smith x[0] -= m1*x1 + m3*x2; 2854eeb42bcSBarry Smith x[1] -= m2*x1 + m4*x2; 2864eeb42bcSBarry Smith x[2] -= m1*x3 + m3*x4; 2874eeb42bcSBarry Smith x[3] -= m2*x3 + m4*x4; 2884eeb42bcSBarry Smith pv += 4; 2894eeb42bcSBarry Smith } 290dc0b31edSSatish Balay ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr); 2914eeb42bcSBarry Smith } 2924eeb42bcSBarry Smith row = *ajtmp++; 2934eeb42bcSBarry Smith } 2944eeb42bcSBarry Smith /* finished row so stick it into b->a */ 2954078e994SBarry Smith pv = ba + 4*bi[i]; 2964078e994SBarry Smith pj = bj + bi[i]; 2974078e994SBarry Smith nz = bi[i+1] - bi[i]; 2984eeb42bcSBarry Smith for (j=0; j<nz; j++) { 2994eeb42bcSBarry Smith x = rtmp+4*pj[j]; 3004eeb42bcSBarry Smith pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 3014eeb42bcSBarry Smith pv += 4; 3024eeb42bcSBarry Smith } 3034eeb42bcSBarry Smith /* invert diagonal block */ 3044078e994SBarry Smith w = ba + 4*diag_offset[i]; 30562bba022SBarry Smith ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr); 3064eeb42bcSBarry Smith } 3074eeb42bcSBarry Smith 308606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 3094eeb42bcSBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 3104eeb42bcSBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 31106e38f1dSHong Zhang C->ops->solve = MatSolve_SeqBAIJ_2_inplace; 31206e38f1dSHong Zhang C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_inplace; 3134eeb42bcSBarry Smith C->assembled = PETSC_TRUE; 314766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 3153a40ed3dSBarry Smith PetscFunctionReturn(0); 3164eeb42bcSBarry Smith } 3174cd38560SBarry Smith /* 3184cd38560SBarry Smith Version for when blocks are 2 by 2 Using natural ordering 3194cd38560SBarry Smith */ 3204a2ae208SSatish Balay #undef __FUNCT__ 32106e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace" 32206e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info) 3234cd38560SBarry Smith { 3244cd38560SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 325dfbe8321SBarry Smith PetscErrorCode ierr; 326690b6cddSBarry Smith PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 327690b6cddSBarry Smith PetscInt *ajtmpold,*ajtmp,nz,row; 328690b6cddSBarry Smith PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; 329329f5518SBarry Smith MatScalar *pv,*v,*rtmp,*pc,*w,*x; 3302515f8d2SSatish Balay MatScalar p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4; 3314cd38560SBarry Smith MatScalar *ba = b->a,*aa = a->a; 332182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 3334cd38560SBarry Smith 3344cd38560SBarry Smith PetscFunctionBegin; 335b0a32e0cSBarry Smith ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 3364cd38560SBarry Smith for (i=0; i<n; i++) { 3374cd38560SBarry Smith nz = bi[i+1] - bi[i]; 3384cd38560SBarry Smith ajtmp = bj + bi[i]; 3394cd38560SBarry Smith for (j=0; j<nz; j++) { 3404cd38560SBarry Smith x = rtmp+4*ajtmp[j]; 3414cd38560SBarry Smith x[0] = x[1] = x[2] = x[3] = 0.0; 3424cd38560SBarry Smith } 3434cd38560SBarry Smith /* load in initial (unfactored row) */ 3444cd38560SBarry Smith nz = ai[i+1] - ai[i]; 3454cd38560SBarry Smith ajtmpold = aj + ai[i]; 3464cd38560SBarry Smith v = aa + 4*ai[i]; 3474cd38560SBarry Smith for (j=0; j<nz; j++) { 3484cd38560SBarry Smith x = rtmp+4*ajtmpold[j]; 3494cd38560SBarry Smith x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 3504cd38560SBarry Smith v += 4; 3514cd38560SBarry Smith } 3524cd38560SBarry Smith row = *ajtmp++; 3534cd38560SBarry Smith while (row < i) { 3544cd38560SBarry Smith pc = rtmp + 4*row; 3554cd38560SBarry Smith p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 3564cd38560SBarry Smith if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) { 3574cd38560SBarry Smith pv = ba + 4*diag_offset[row]; 3584cd38560SBarry Smith pj = bj + diag_offset[row] + 1; 3594cd38560SBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 3604cd38560SBarry Smith pc[0] = m1 = p1*x1 + p3*x2; 3614cd38560SBarry Smith pc[1] = m2 = p2*x1 + p4*x2; 3624cd38560SBarry Smith pc[2] = m3 = p1*x3 + p3*x4; 3634cd38560SBarry Smith pc[3] = m4 = p2*x3 + p4*x4; 3644cd38560SBarry Smith nz = bi[row+1] - diag_offset[row] - 1; 3654cd38560SBarry Smith pv += 4; 3664cd38560SBarry Smith for (j=0; j<nz; j++) { 3674cd38560SBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 3684cd38560SBarry Smith x = rtmp + 4*pj[j]; 3694cd38560SBarry Smith x[0] -= m1*x1 + m3*x2; 3704cd38560SBarry Smith x[1] -= m2*x1 + m4*x2; 3714cd38560SBarry Smith x[2] -= m1*x3 + m3*x4; 3724cd38560SBarry Smith x[3] -= m2*x3 + m4*x4; 3734cd38560SBarry Smith pv += 4; 3744cd38560SBarry Smith } 375dc0b31edSSatish Balay ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr); 3764cd38560SBarry Smith } 3774cd38560SBarry Smith row = *ajtmp++; 3784cd38560SBarry Smith } 3794cd38560SBarry Smith /* finished row so stick it into b->a */ 3804cd38560SBarry Smith pv = ba + 4*bi[i]; 3814cd38560SBarry Smith pj = bj + bi[i]; 3824cd38560SBarry Smith nz = bi[i+1] - bi[i]; 3834cd38560SBarry Smith for (j=0; j<nz; j++) { 3844cd38560SBarry Smith x = rtmp+4*pj[j]; 3854cd38560SBarry Smith pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 3862efa7f71SHong Zhang /* 3872efa7f71SHong Zhang printf(" col %d:",pj[j]); 3882efa7f71SHong Zhang PetscInt j1; 3892efa7f71SHong Zhang for (j1=0; j1<4; j1++) printf(" %g,",*(pv+j1)); 3902efa7f71SHong Zhang printf("\n"); 3912efa7f71SHong Zhang */ 3924cd38560SBarry Smith pv += 4; 3934cd38560SBarry Smith } 3944cd38560SBarry Smith /* invert diagonal block */ 3954cd38560SBarry Smith w = ba + 4*diag_offset[i]; 3962efa7f71SHong Zhang /* 3972efa7f71SHong Zhang printf(" \n%d -th: diag: ",i); 3982efa7f71SHong Zhang for (j=0; j<4; j++){ 3992efa7f71SHong Zhang printf(" %g,",w[j]); 4002efa7f71SHong Zhang } 4012efa7f71SHong Zhang printf("\n----------------------------\n"); 4022efa7f71SHong Zhang */ 40362bba022SBarry Smith ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr); 4044cd38560SBarry Smith } 4054cd38560SBarry Smith 406606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 40706e38f1dSHong Zhang C->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering_inplace; 40806e38f1dSHong Zhang C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering_inplace; 4094cd38560SBarry Smith C->assembled = PETSC_TRUE; 410766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 4114cd38560SBarry Smith PetscFunctionReturn(0); 4124cd38560SBarry Smith } 4137fc0212eSBarry Smith 4147fc0212eSBarry Smith /* ----------------------------------------------------------- */ 4157fc0212eSBarry Smith /* 4167fc0212eSBarry Smith Version for when blocks are 1 by 1. 4177fc0212eSBarry Smith */ 4184a2ae208SSatish Balay #undef __FUNCT__ 419*048b5e81SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1" 420*048b5e81SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat B,Mat A,const MatFactorInfo *info) 421*048b5e81SShri Abhyankar { 422*048b5e81SShri Abhyankar Mat C=B; 423*048b5e81SShri Abhyankar Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; 424*048b5e81SShri Abhyankar IS isrow = b->row,isicol = b->icol; 425*048b5e81SShri Abhyankar PetscErrorCode ierr; 426*048b5e81SShri Abhyankar const PetscInt *r,*ic,*ics; 427*048b5e81SShri Abhyankar const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bdiag=b->diag; 428*048b5e81SShri Abhyankar PetscInt i,j,k,nz,nzL,row,*pj; 429*048b5e81SShri Abhyankar const PetscInt *ajtmp,*bjtmp; 430*048b5e81SShri Abhyankar MatScalar *rtmp,*pc,multiplier,*pv; 431*048b5e81SShri Abhyankar const MatScalar *aa=a->a,*v; 432*048b5e81SShri Abhyankar PetscTruth row_identity,col_identity; 433*048b5e81SShri Abhyankar FactorShiftCtx sctx; 434*048b5e81SShri Abhyankar const PetscInt *ddiag; 435*048b5e81SShri Abhyankar PetscReal rs; 436*048b5e81SShri Abhyankar MatScalar d; 437*048b5e81SShri Abhyankar 438*048b5e81SShri Abhyankar PetscFunctionBegin; 439*048b5e81SShri Abhyankar /* MatPivotSetUp(): initialize shift context sctx */ 440*048b5e81SShri Abhyankar ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr); 441*048b5e81SShri Abhyankar 442*048b5e81SShri Abhyankar if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */ 443*048b5e81SShri Abhyankar ddiag = a->diag; 444*048b5e81SShri Abhyankar sctx.shift_top = info->zeropivot; 445*048b5e81SShri Abhyankar for (i=0; i<n; i++) { 446*048b5e81SShri Abhyankar /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */ 447*048b5e81SShri Abhyankar d = (aa)[ddiag[i]]; 448*048b5e81SShri Abhyankar rs = -PetscAbsScalar(d) - PetscRealPart(d); 449*048b5e81SShri Abhyankar v = aa+ai[i]; 450*048b5e81SShri Abhyankar nz = ai[i+1] - ai[i]; 451*048b5e81SShri Abhyankar for (j=0; j<nz; j++) 452*048b5e81SShri Abhyankar rs += PetscAbsScalar(v[j]); 453*048b5e81SShri Abhyankar if (rs>sctx.shift_top) sctx.shift_top = rs; 454*048b5e81SShri Abhyankar } 455*048b5e81SShri Abhyankar sctx.shift_top *= 1.1; 456*048b5e81SShri Abhyankar sctx.nshift_max = 5; 457*048b5e81SShri Abhyankar sctx.shift_lo = 0.; 458*048b5e81SShri Abhyankar sctx.shift_hi = 1.; 459*048b5e81SShri Abhyankar } 460*048b5e81SShri Abhyankar 461*048b5e81SShri Abhyankar ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 462*048b5e81SShri Abhyankar ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 463*048b5e81SShri Abhyankar ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 464*048b5e81SShri Abhyankar ics = ic; 465*048b5e81SShri Abhyankar 466*048b5e81SShri Abhyankar do { 467*048b5e81SShri Abhyankar sctx.newshift = PETSC_FALSE; 468*048b5e81SShri Abhyankar for (i=0; i<n; i++){ 469*048b5e81SShri Abhyankar /* zero rtmp */ 470*048b5e81SShri Abhyankar /* L part */ 471*048b5e81SShri Abhyankar nz = bi[i+1] - bi[i]; 472*048b5e81SShri Abhyankar bjtmp = bj + bi[i]; 473*048b5e81SShri Abhyankar for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0; 474*048b5e81SShri Abhyankar 475*048b5e81SShri Abhyankar /* U part */ 476*048b5e81SShri Abhyankar nz = bdiag[i]-bdiag[i+1]; 477*048b5e81SShri Abhyankar bjtmp = bj + bdiag[i+1]+1; 478*048b5e81SShri Abhyankar for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0; 479*048b5e81SShri Abhyankar 480*048b5e81SShri Abhyankar /* load in initial (unfactored row) */ 481*048b5e81SShri Abhyankar nz = ai[r[i]+1] - ai[r[i]]; 482*048b5e81SShri Abhyankar ajtmp = aj + ai[r[i]]; 483*048b5e81SShri Abhyankar v = aa + ai[r[i]]; 484*048b5e81SShri Abhyankar for (j=0; j<nz; j++) { 485*048b5e81SShri Abhyankar rtmp[ics[ajtmp[j]]] = v[j]; 486*048b5e81SShri Abhyankar } 487*048b5e81SShri Abhyankar /* ZeropivotApply() */ 488*048b5e81SShri Abhyankar rtmp[i] += sctx.shift_amount; /* shift the diagonal of the matrix */ 489*048b5e81SShri Abhyankar 490*048b5e81SShri Abhyankar /* elimination */ 491*048b5e81SShri Abhyankar bjtmp = bj + bi[i]; 492*048b5e81SShri Abhyankar row = *bjtmp++; 493*048b5e81SShri Abhyankar nzL = bi[i+1] - bi[i]; 494*048b5e81SShri Abhyankar for(k=0; k < nzL;k++) { 495*048b5e81SShri Abhyankar pc = rtmp + row; 496*048b5e81SShri Abhyankar if (*pc != 0.0) { 497*048b5e81SShri Abhyankar pv = b->a + bdiag[row]; 498*048b5e81SShri Abhyankar multiplier = *pc * (*pv); 499*048b5e81SShri Abhyankar *pc = multiplier; 500*048b5e81SShri Abhyankar pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */ 501*048b5e81SShri Abhyankar pv = b->a + bdiag[row+1]+1; 502*048b5e81SShri Abhyankar nz = bdiag[row]-bdiag[row+1]-1; /* num of entries in U(row,:) excluding diag */ 503*048b5e81SShri Abhyankar for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j]; 504*048b5e81SShri Abhyankar ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); 505*048b5e81SShri Abhyankar } 506*048b5e81SShri Abhyankar row = *bjtmp++; 507*048b5e81SShri Abhyankar } 508*048b5e81SShri Abhyankar 509*048b5e81SShri Abhyankar /* finished row so stick it into b->a */ 510*048b5e81SShri Abhyankar rs = 0.0; 511*048b5e81SShri Abhyankar /* L part */ 512*048b5e81SShri Abhyankar pv = b->a + bi[i] ; 513*048b5e81SShri Abhyankar pj = b->j + bi[i] ; 514*048b5e81SShri Abhyankar nz = bi[i+1] - bi[i]; 515*048b5e81SShri Abhyankar for (j=0; j<nz; j++) { 516*048b5e81SShri Abhyankar pv[j] = rtmp[pj[j]]; rs += PetscAbsScalar(pv[j]); 517*048b5e81SShri Abhyankar } 518*048b5e81SShri Abhyankar 519*048b5e81SShri Abhyankar /* U part */ 520*048b5e81SShri Abhyankar pv = b->a + bdiag[i+1]+1; 521*048b5e81SShri Abhyankar pj = b->j + bdiag[i+1]+1; 522*048b5e81SShri Abhyankar nz = bdiag[i] - bdiag[i+1]-1; 523*048b5e81SShri Abhyankar for (j=0; j<nz; j++) { 524*048b5e81SShri Abhyankar pv[j] = rtmp[pj[j]]; rs += PetscAbsScalar(pv[j]); 525*048b5e81SShri Abhyankar } 526*048b5e81SShri Abhyankar 527*048b5e81SShri Abhyankar sctx.rs = rs; 528*048b5e81SShri Abhyankar sctx.pv = rtmp[i]; 529*048b5e81SShri Abhyankar ierr = MatPivotCheck(info,&sctx,i);CHKERRQ(ierr); 530*048b5e81SShri Abhyankar if(sctx.newshift) break; /* break for-loop */ 531*048b5e81SShri Abhyankar rtmp[i] = sctx.pv; /* sctx.pv might be updated in the case of MAT_SHIFT_INBLOCKS */ 532*048b5e81SShri Abhyankar 533*048b5e81SShri Abhyankar /* Mark diagonal and invert diagonal for simplier triangular solves */ 534*048b5e81SShri Abhyankar pv = b->a + bdiag[i]; 535*048b5e81SShri Abhyankar *pv = 1.0/rtmp[i]; 536*048b5e81SShri Abhyankar 537*048b5e81SShri Abhyankar } /* endof for (i=0; i<n; i++){ */ 538*048b5e81SShri Abhyankar 539*048b5e81SShri Abhyankar /* MatPivotRefine() */ 540*048b5e81SShri Abhyankar if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE && !sctx.newshift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max){ 541*048b5e81SShri Abhyankar /* 542*048b5e81SShri Abhyankar * if no shift in this attempt & shifting & started shifting & can refine, 543*048b5e81SShri Abhyankar * then try lower shift 544*048b5e81SShri Abhyankar */ 545*048b5e81SShri Abhyankar sctx.shift_hi = sctx.shift_fraction; 546*048b5e81SShri Abhyankar sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.; 547*048b5e81SShri Abhyankar sctx.shift_amount = sctx.shift_fraction * sctx.shift_top; 548*048b5e81SShri Abhyankar sctx.newshift = PETSC_TRUE; 549*048b5e81SShri Abhyankar sctx.nshift++; 550*048b5e81SShri Abhyankar } 551*048b5e81SShri Abhyankar } while (sctx.newshift); 552*048b5e81SShri Abhyankar 553*048b5e81SShri Abhyankar ierr = PetscFree(rtmp);CHKERRQ(ierr); 554*048b5e81SShri Abhyankar ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 555*048b5e81SShri Abhyankar ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 556*048b5e81SShri Abhyankar 557*048b5e81SShri Abhyankar ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 558*048b5e81SShri Abhyankar ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr); 559*048b5e81SShri Abhyankar if (row_identity && col_identity) { 560*048b5e81SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_1_NaturalOrdering; 561*048b5e81SShri Abhyankar } else { 562*048b5e81SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_1; 563*048b5e81SShri Abhyankar } 564*048b5e81SShri Abhyankar C->ops->solveadd = MatSolveAdd_SeqAIJ; 565*048b5e81SShri Abhyankar C->assembled = PETSC_TRUE; 566*048b5e81SShri Abhyankar C->preallocated = PETSC_TRUE; 567*048b5e81SShri Abhyankar ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr); 568*048b5e81SShri Abhyankar 569*048b5e81SShri Abhyankar /* MatShiftView(A,info,&sctx) */ 570*048b5e81SShri Abhyankar if (sctx.nshift){ 571*048b5e81SShri Abhyankar if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { 572*048b5e81SShri Abhyankar ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,sctx.shift_fraction,sctx.shift_top);CHKERRQ(ierr); 573*048b5e81SShri Abhyankar } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 574*048b5e81SShri Abhyankar ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 575*048b5e81SShri Abhyankar } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS){ 576*048b5e81SShri Abhyankar ierr = PetscInfo2(A,"number of shift_inblocks applied %D, each shift_amount %G\n",sctx.nshift,info->shiftamount);CHKERRQ(ierr); 577*048b5e81SShri Abhyankar } 578*048b5e81SShri Abhyankar } 579*048b5e81SShri Abhyankar PetscFunctionReturn(0); 580*048b5e81SShri Abhyankar } 581*048b5e81SShri Abhyankar 582*048b5e81SShri Abhyankar #undef __FUNCT__ 58306e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1_inplace" 58406e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info) 5857fc0212eSBarry Smith { 5867fc0212eSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 5877cdf2dbbSSatish Balay IS isrow = b->row,isicol = b->icol; 5886849ba73SBarry Smith PetscErrorCode ierr; 5895d0c19d7SBarry Smith const PetscInt *r,*ic; 5905d0c19d7SBarry Smith PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 591690b6cddSBarry Smith PetscInt *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j; 592690b6cddSBarry Smith PetscInt *diag_offset = b->diag,diag,*pj; 593329f5518SBarry Smith MatScalar *pv,*v,*rtmp,multiplier,*pc; 5943f1db9ecSBarry Smith MatScalar *ba = b->a,*aa = a->a; 595f58c8c74SBarry Smith PetscTruth row_identity, col_identity; 5967fc0212eSBarry Smith 5973a40ed3dSBarry Smith PetscFunctionBegin; 5987fc0212eSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 5997fc0212eSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 600b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 6017fc0212eSBarry Smith 6027fc0212eSBarry Smith for (i=0; i<n; i++) { 6034078e994SBarry Smith nz = bi[i+1] - bi[i]; 6044078e994SBarry Smith ajtmp = bj + bi[i]; 6057fc0212eSBarry Smith for (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0; 6067fc0212eSBarry Smith 6077fc0212eSBarry Smith /* load in initial (unfactored row) */ 6084078e994SBarry Smith nz = ai[r[i]+1] - ai[r[i]]; 6094078e994SBarry Smith ajtmpold = aj + ai[r[i]]; 6104078e994SBarry Smith v = aa + ai[r[i]]; 6117fc0212eSBarry Smith for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] = v[j]; 6127fc0212eSBarry Smith 6137fc0212eSBarry Smith row = *ajtmp++; 6147fc0212eSBarry Smith while (row < i) { 6157fc0212eSBarry Smith pc = rtmp + row; 6167fc0212eSBarry Smith if (*pc != 0.0) { 6174078e994SBarry Smith pv = ba + diag_offset[row]; 6184078e994SBarry Smith pj = bj + diag_offset[row] + 1; 6197fc0212eSBarry Smith multiplier = *pc * *pv++; 6207fc0212eSBarry Smith *pc = multiplier; 6214078e994SBarry Smith nz = bi[row+1] - diag_offset[row] - 1; 6227fc0212eSBarry Smith for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j]; 623dc0b31edSSatish Balay ierr = PetscLogFlops(1.0+2.0*nz);CHKERRQ(ierr); 6247fc0212eSBarry Smith } 6257fc0212eSBarry Smith row = *ajtmp++; 6267fc0212eSBarry Smith } 6277fc0212eSBarry Smith /* finished row so stick it into b->a */ 6284078e994SBarry Smith pv = ba + bi[i]; 6294078e994SBarry Smith pj = bj + bi[i]; 6304078e994SBarry Smith nz = bi[i+1] - bi[i]; 6317fc0212eSBarry Smith for (j=0; j<nz; j++) {pv[j] = rtmp[pj[j]];} 6324078e994SBarry Smith diag = diag_offset[i] - bi[i]; 6337fc0212eSBarry Smith /* check pivot entry for current row */ 63465e19b50SBarry Smith if (pv[diag] == 0.0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot: row in original ordering %D in permuted ordering %D",r[i],i); 6357fc0212eSBarry Smith pv[diag] = 1.0/pv[diag]; 6367fc0212eSBarry Smith } 6377fc0212eSBarry Smith 638606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 6397fc0212eSBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 6407fc0212eSBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 641f58c8c74SBarry Smith ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 642f58c8c74SBarry Smith ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr); 643f58c8c74SBarry Smith if (row_identity && col_identity) { 64406e38f1dSHong Zhang C->ops->solve = MatSolve_SeqBAIJ_1_NaturalOrdering_inplace; 64506e38f1dSHong Zhang C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering_inplace; 646f58c8c74SBarry Smith } else { 64706e38f1dSHong Zhang C->ops->solve = MatSolve_SeqBAIJ_1_inplace; 64806e38f1dSHong Zhang C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_inplace; 649f58c8c74SBarry Smith } 6507fc0212eSBarry Smith C->assembled = PETSC_TRUE; 651d0f46423SBarry Smith ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr); 6523a40ed3dSBarry Smith PetscFunctionReturn(0); 6537fc0212eSBarry Smith } 6547fc0212eSBarry Smith 655e631078cSBarry Smith EXTERN_C_BEGIN 656b24902e0SBarry Smith #undef __FUNCT__ 657b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqbaij_petsc" 658b24902e0SBarry Smith PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A,MatFactorType ftype,Mat *B) 659b24902e0SBarry Smith { 660d0f46423SBarry Smith PetscInt n = A->rmap->n; 661b24902e0SBarry Smith PetscErrorCode ierr; 662b24902e0SBarry Smith 663b24902e0SBarry Smith PetscFunctionBegin; 664b24902e0SBarry Smith ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr); 665b24902e0SBarry Smith ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 666c8342467SHong Zhang if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) { 667b24902e0SBarry Smith ierr = MatSetType(*B,MATSEQBAIJ);CHKERRQ(ierr); 6684dd39f65SShri Abhyankar (*B)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqBAIJ; 6694dd39f65SShri Abhyankar (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqBAIJ; 670b24902e0SBarry Smith } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 671b24902e0SBarry Smith ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr); 6725c9eb25fSBarry Smith ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 6735c9eb25fSBarry Smith (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqBAIJ; 6745c9eb25fSBarry Smith (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ; 675e32f2f54SBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 676d5f3da31SBarry Smith (*B)->factortype = ftype; 677b24902e0SBarry Smith PetscFunctionReturn(0); 678b24902e0SBarry Smith } 679e631078cSBarry Smith EXTERN_C_END 680273d9f13SBarry Smith 681db4efbfdSBarry Smith EXTERN_C_BEGIN 682db4efbfdSBarry Smith #undef __FUNCT__ 683db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactorAvailable_seqaij_petsc" 684db4efbfdSBarry Smith PetscErrorCode MatGetFactorAvailable_seqbaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg) 685db4efbfdSBarry Smith { 686db4efbfdSBarry Smith PetscFunctionBegin; 687db4efbfdSBarry Smith *flg = PETSC_TRUE; 688db4efbfdSBarry Smith PetscFunctionReturn(0); 689db4efbfdSBarry Smith } 690db4efbfdSBarry Smith EXTERN_C_END 691db4efbfdSBarry Smith 6925d517e7eSBarry Smith /* ----------------------------------------------------------- */ 6934a2ae208SSatish Balay #undef __FUNCT__ 6944a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqBAIJ" 6950481f469SBarry Smith PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,const MatFactorInfo *info) 6965d517e7eSBarry Smith { 697dfbe8321SBarry Smith PetscErrorCode ierr; 6985d517e7eSBarry Smith Mat C; 6995d517e7eSBarry Smith 7003a40ed3dSBarry Smith PetscFunctionBegin; 70143244d56SBarry Smith ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr); 702719d5645SBarry Smith ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr); 703719d5645SBarry Smith ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr); 704db4efbfdSBarry Smith A->ops->solve = C->ops->solve; 705db4efbfdSBarry Smith A->ops->solvetranspose = C->ops->solvetranspose; 706273d9f13SBarry Smith ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 70752e6d16bSBarry Smith ierr = PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr); 7083a40ed3dSBarry Smith PetscFunctionReturn(0); 7095d517e7eSBarry Smith } 7104cd38560SBarry Smith 7117c4f633dSBarry Smith #include "../src/mat/impls/sbaij/seq/sbaij.h" 712c05c3958SHong Zhang #undef __FUNCT__ 713c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N" 7140481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat C,Mat A,const MatFactorInfo *info) 715c05c3958SHong Zhang { 716f3086b4bSHong Zhang PetscErrorCode ierr; 71778910aadSHong Zhang Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 71878910aadSHong Zhang Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 71978910aadSHong Zhang IS ip=b->row; 7205d0c19d7SBarry Smith const PetscInt *rip; 7215d0c19d7SBarry Smith PetscInt i,j,mbs=a->mbs,bs=A->rmap->bs,*bi=b->i,*bj=b->j,*bcol; 72278910aadSHong Zhang PetscInt *ai=a->i,*aj=a->j; 72378910aadSHong Zhang PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 72478910aadSHong Zhang MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 72507b50cabSHong Zhang PetscReal rs; 7260e95ead3SHong Zhang FactorShiftCtx sctx; 72778910aadSHong Zhang 728c05c3958SHong Zhang PetscFunctionBegin; 72907b50cabSHong Zhang if (bs > 1) {/* convert A to a SBAIJ matrix and apply Cholesky factorization from it */ 7306ad2eaddSHong Zhang if (!a->sbaijMat){ 731ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 7326ad2eaddSHong Zhang } 733719d5645SBarry Smith ierr = (a->sbaijMat)->ops->choleskyfactornumeric(C,a->sbaijMat,info);CHKERRQ(ierr); 7346ad2eaddSHong Zhang ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr); 7356ad2eaddSHong Zhang a->sbaijMat = PETSC_NULL; 7366ad2eaddSHong Zhang PetscFunctionReturn(0); 7376ad2eaddSHong Zhang } 73878910aadSHong Zhang 73907b50cabSHong Zhang /* MatPivotSetUp(): initialize shift context sctx */ 74007b50cabSHong Zhang ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr); 7413cea4cbeSHong Zhang 7426ad2eaddSHong Zhang ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 743fca92195SBarry Smith ierr = PetscMalloc3(mbs,MatScalar,&rtmp,mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr); 74478910aadSHong Zhang 74575567043SBarry Smith sctx.shift_amount = 0.; 7463cea4cbeSHong Zhang sctx.nshift = 0; 74778910aadSHong Zhang do { 74807b50cabSHong Zhang sctx.newshift = PETSC_FALSE; 74978910aadSHong Zhang for (i=0; i<mbs; i++) { 750e60cf9a0SBarry Smith rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 75178910aadSHong Zhang } 75278910aadSHong Zhang 75378910aadSHong Zhang for (k = 0; k<mbs; k++){ 75478910aadSHong Zhang bval = ba + bi[k]; 75578910aadSHong Zhang /* initialize k-th row by the perm[k]-th row of A */ 75678910aadSHong Zhang jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 75778910aadSHong Zhang for (j = jmin; j < jmax; j++){ 75878910aadSHong Zhang col = rip[aj[j]]; 75978910aadSHong Zhang if (col >= k){ /* only take upper triangular entry */ 76078910aadSHong Zhang rtmp[col] = aa[j]; 76178910aadSHong Zhang *bval++ = 0.0; /* for in-place factorization */ 76278910aadSHong Zhang } 76378910aadSHong Zhang } 7643cea4cbeSHong Zhang 7653cea4cbeSHong Zhang /* shift the diagonal of the matrix */ 7663cea4cbeSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 76778910aadSHong Zhang 76878910aadSHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 76978910aadSHong Zhang dk = rtmp[k]; 77078910aadSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 77178910aadSHong Zhang 77278910aadSHong Zhang while (i < k){ 77378910aadSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 77478910aadSHong Zhang 77578910aadSHong Zhang /* compute multiplier, update diag(k) and U(i,k) */ 77678910aadSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 77778910aadSHong Zhang uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 77878910aadSHong Zhang dk += uikdi*ba[ili]; 77978910aadSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 78078910aadSHong Zhang 78178910aadSHong Zhang /* add multiple of row i to k-th row */ 78278910aadSHong Zhang jmin = ili + 1; jmax = bi[i+1]; 78378910aadSHong Zhang if (jmin < jmax){ 78478910aadSHong Zhang for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 78578910aadSHong Zhang /* update il and jl for row i */ 78678910aadSHong Zhang il[i] = jmin; 78778910aadSHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 78878910aadSHong Zhang } 78978910aadSHong Zhang i = nexti; 79078910aadSHong Zhang } 79178910aadSHong Zhang 7923cea4cbeSHong Zhang /* shift the diagonals when zero pivot is detected */ 7933cea4cbeSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 7943cea4cbeSHong Zhang rs = 0.0; 79578910aadSHong Zhang jmin = bi[k]+1; 79678910aadSHong Zhang nz = bi[k+1] - jmin; 79778910aadSHong Zhang if (nz){ 79878910aadSHong Zhang bcol = bj + jmin; 79978910aadSHong Zhang while (nz--){ 80089f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 80189f845c8SHong Zhang bcol++; 80278910aadSHong Zhang } 80378910aadSHong Zhang } 8043cea4cbeSHong Zhang 8053cea4cbeSHong Zhang sctx.rs = rs; 8063cea4cbeSHong Zhang sctx.pv = dk; 80707b50cabSHong Zhang ierr = MatPivotCheck(info,&sctx,k);CHKERRQ(ierr); 80807b50cabSHong Zhang if (sctx.newshift) break; 8090e95ead3SHong Zhang dk = sctx.pv; 81078910aadSHong Zhang 81178910aadSHong Zhang /* copy data into U(k,:) */ 81278910aadSHong Zhang ba[bi[k]] = 1.0/dk; /* U(k,k) */ 81378910aadSHong Zhang jmin = bi[k]+1; jmax = bi[k+1]; 81478910aadSHong Zhang if (jmin < jmax) { 81578910aadSHong Zhang for (j=jmin; j<jmax; j++){ 81678910aadSHong Zhang col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 81778910aadSHong Zhang } 81878910aadSHong Zhang /* add the k-th row into il and jl */ 81978910aadSHong Zhang il[k] = jmin; 82078910aadSHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 82178910aadSHong Zhang } 82278910aadSHong Zhang } 82307b50cabSHong Zhang } while (sctx.newshift); 824fca92195SBarry Smith ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr); 82578910aadSHong Zhang 82678910aadSHong Zhang ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 82778910aadSHong Zhang C->assembled = PETSC_TRUE; 82878910aadSHong Zhang C->preallocated = PETSC_TRUE; 829d0f46423SBarry Smith ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); 8303cea4cbeSHong Zhang if (sctx.nshift){ 831f4db908eSBarry Smith if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { 8321e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 833f4db908eSBarry Smith } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 834e738e422SBarry Smith ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 83578910aadSHong Zhang } 83678910aadSHong Zhang } 837c05c3958SHong Zhang PetscFunctionReturn(0); 838c05c3958SHong Zhang } 8394cd38560SBarry Smith 840c05c3958SHong Zhang #undef __FUNCT__ 841c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering" 8420481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info) 843c05c3958SHong Zhang { 84478910aadSHong Zhang Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 84578910aadSHong Zhang Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 84678910aadSHong Zhang PetscErrorCode ierr; 8473cea4cbeSHong Zhang PetscInt i,j,am=a->mbs; 84878910aadSHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 8493cea4cbeSHong Zhang PetscInt k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz; 85078910aadSHong Zhang MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval; 8510e95ead3SHong Zhang PetscReal rs; 8520e95ead3SHong Zhang FactorShiftCtx sctx; 85378910aadSHong Zhang 854c05c3958SHong Zhang PetscFunctionBegin; 8550e95ead3SHong Zhang /* MatPivotSetUp(): initialize shift context sctx */ 8560e95ead3SHong Zhang ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr); 8573cea4cbeSHong Zhang 858fca92195SBarry Smith ierr = PetscMalloc3(am,MatScalar,&rtmp,am,PetscInt,&il,am,PetscInt,&jl);CHKERRQ(ierr); 85978910aadSHong Zhang 86078910aadSHong Zhang do { 86107b50cabSHong Zhang sctx.newshift = PETSC_FALSE; 86278910aadSHong Zhang for (i=0; i<am; i++) { 863e60cf9a0SBarry Smith rtmp[i] = 0.0; jl[i] = am; il[0] = 0; 86478910aadSHong Zhang } 86578910aadSHong Zhang 86678910aadSHong Zhang for (k = 0; k<am; k++){ 86778910aadSHong Zhang /* initialize k-th row with elements nonzero in row perm(k) of A */ 86878910aadSHong Zhang nz = ai[k+1] - ai[k]; 86978910aadSHong Zhang acol = aj + ai[k]; 87078910aadSHong Zhang aval = aa + ai[k]; 87178910aadSHong Zhang bval = ba + bi[k]; 87278910aadSHong Zhang while (nz -- ){ 87378910aadSHong Zhang if (*acol < k) { /* skip lower triangular entries */ 87478910aadSHong Zhang acol++; aval++; 87578910aadSHong Zhang } else { 87678910aadSHong Zhang rtmp[*acol++] = *aval++; 87778910aadSHong Zhang *bval++ = 0.0; /* for in-place factorization */ 87878910aadSHong Zhang } 87978910aadSHong Zhang } 8803cea4cbeSHong Zhang 8813cea4cbeSHong Zhang /* shift the diagonal of the matrix */ 8823cea4cbeSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 88378910aadSHong Zhang 88478910aadSHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 88578910aadSHong Zhang dk = rtmp[k]; 88678910aadSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 88778910aadSHong Zhang 88878910aadSHong Zhang while (i < k){ 88978910aadSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 89078910aadSHong Zhang /* compute multiplier, update D(k) and U(i,k) */ 89178910aadSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 89278910aadSHong Zhang uikdi = - ba[ili]*ba[bi[i]]; 89378910aadSHong Zhang dk += uikdi*ba[ili]; 89478910aadSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 89578910aadSHong Zhang 89678910aadSHong Zhang /* add multiple of row i to k-th row ... */ 89778910aadSHong Zhang jmin = ili + 1; 89878910aadSHong Zhang nz = bi[i+1] - jmin; 89978910aadSHong Zhang if (nz > 0){ 90078910aadSHong Zhang bcol = bj + jmin; 90178910aadSHong Zhang bval = ba + jmin; 90278910aadSHong Zhang while (nz --) rtmp[*bcol++] += uikdi*(*bval++); 90378910aadSHong Zhang /* update il and jl for i-th row */ 90478910aadSHong Zhang il[i] = jmin; 90578910aadSHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 90678910aadSHong Zhang } 90778910aadSHong Zhang i = nexti; 90878910aadSHong Zhang } 90978910aadSHong Zhang 9103cea4cbeSHong Zhang /* shift the diagonals when zero pivot is detected */ 9113cea4cbeSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 9123cea4cbeSHong Zhang rs = 0.0; 91378910aadSHong Zhang jmin = bi[k]+1; 91478910aadSHong Zhang nz = bi[k+1] - jmin; 91578910aadSHong Zhang if (nz){ 91678910aadSHong Zhang bcol = bj + jmin; 91778910aadSHong Zhang while (nz--){ 91889f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 91989f845c8SHong Zhang bcol++; 92078910aadSHong Zhang } 92178910aadSHong Zhang } 9223cea4cbeSHong Zhang 9233cea4cbeSHong Zhang sctx.rs = rs; 9243cea4cbeSHong Zhang sctx.pv = dk; 92507b50cabSHong Zhang ierr = MatPivotCheck(info,&sctx,k);CHKERRQ(ierr); 92607b50cabSHong Zhang if (sctx.newshift) break; /* sctx.shift_amount is updated */ 9270e95ead3SHong Zhang dk = sctx.pv; 92878910aadSHong Zhang 92978910aadSHong Zhang /* copy data into U(k,:) */ 93078910aadSHong Zhang ba[bi[k]] = 1.0/dk; 93178910aadSHong Zhang jmin = bi[k]+1; 93278910aadSHong Zhang nz = bi[k+1] - jmin; 93378910aadSHong Zhang if (nz){ 93478910aadSHong Zhang bcol = bj + jmin; 93578910aadSHong Zhang bval = ba + jmin; 93678910aadSHong Zhang while (nz--){ 93778910aadSHong Zhang *bval++ = rtmp[*bcol]; 93878910aadSHong Zhang rtmp[*bcol++] = 0.0; 93978910aadSHong Zhang } 94078910aadSHong Zhang /* add k-th row into il and jl */ 94178910aadSHong Zhang il[k] = jmin; 94278910aadSHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 94378910aadSHong Zhang } 94478910aadSHong Zhang } 94507b50cabSHong Zhang } while (sctx.newshift); 946fca92195SBarry Smith ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr); 94778910aadSHong Zhang 9480a3c351aSHong Zhang C->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 9490a3c351aSHong Zhang C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 95078910aadSHong Zhang C->assembled = PETSC_TRUE; 95178910aadSHong Zhang C->preallocated = PETSC_TRUE; 952d0f46423SBarry Smith ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); 9533cea4cbeSHong Zhang if (sctx.nshift){ 954f4db908eSBarry Smith if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 9551e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 956f4db908eSBarry Smith } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { 9571e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 95878910aadSHong Zhang } 95978910aadSHong Zhang } 960c05c3958SHong Zhang PetscFunctionReturn(0); 961c05c3958SHong Zhang } 962c05c3958SHong Zhang 963c05c3958SHong Zhang #include "petscbt.h" 9647c4f633dSBarry Smith #include "../src/mat/utils/freespace.h" 965c05c3958SHong Zhang #undef __FUNCT__ 966c05c3958SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ" 9670481f469SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 968c05c3958SHong Zhang { 96978910aadSHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 97078910aadSHong Zhang Mat_SeqSBAIJ *b; 97178910aadSHong Zhang Mat B; 97278910aadSHong Zhang PetscErrorCode ierr; 97378910aadSHong Zhang PetscTruth perm_identity; 9745d0c19d7SBarry Smith PetscInt reallocs=0,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap->bs,*ui; 9755d0c19d7SBarry Smith const PetscInt *rip; 97678910aadSHong Zhang PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 977cfdb8b8aSHong Zhang PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr; 97878910aadSHong Zhang PetscReal fill=info->fill,levels=info->levels; 979a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 980a1a86e44SBarry Smith PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 98178910aadSHong Zhang PetscBT lnkbt; 98278910aadSHong Zhang 983c05c3958SHong Zhang PetscFunctionBegin; 9846ad2eaddSHong Zhang if (bs > 1){ 9856ad2eaddSHong Zhang if (!a->sbaijMat){ 986ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 9876ad2eaddSHong Zhang } 988719d5645SBarry Smith (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */ 989719d5645SBarry Smith ierr = MatICCFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr); 9906ad2eaddSHong Zhang PetscFunctionReturn(0); 9916ad2eaddSHong Zhang } 9926ad2eaddSHong Zhang 99378910aadSHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 99478910aadSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 99578910aadSHong Zhang 99678910aadSHong Zhang /* special case that simply copies fill pattern */ 99778910aadSHong Zhang if (!levels && perm_identity) { 99878910aadSHong Zhang ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr); 99978910aadSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 100078910aadSHong Zhang for (i=0; i<am; i++) { 100178910aadSHong Zhang ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */ 100278910aadSHong Zhang } 1003719d5645SBarry Smith B = fact; 100478910aadSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr); 100578910aadSHong Zhang 1006b24902e0SBarry Smith 100778910aadSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 100878910aadSHong Zhang uj = b->j; 100978910aadSHong Zhang for (i=0; i<am; i++) { 101078910aadSHong Zhang aj = a->j + a->diag[i]; 101178910aadSHong Zhang for (j=0; j<ui[i]; j++){ 101278910aadSHong Zhang *uj++ = *aj++; 101378910aadSHong Zhang } 101478910aadSHong Zhang b->ilen[i] = ui[i]; 101578910aadSHong Zhang } 101678910aadSHong Zhang ierr = PetscFree(ui);CHKERRQ(ierr); 1017d5f3da31SBarry Smith B->factortype = MAT_FACTOR_NONE; 101878910aadSHong Zhang ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 101978910aadSHong Zhang ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1020d5f3da31SBarry Smith B->factortype = MAT_FACTOR_ICC; 102178910aadSHong Zhang 102278910aadSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 102378910aadSHong Zhang PetscFunctionReturn(0); 102478910aadSHong Zhang } 102578910aadSHong Zhang 102678910aadSHong Zhang /* initialization */ 102778910aadSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 1028e60cf9a0SBarry Smith ui[0] = 0; 102978910aadSHong Zhang ierr = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr); 103078910aadSHong Zhang 103178910aadSHong Zhang /* jl: linked list for storing indices of the pivot rows 103278910aadSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 1033fca92195SBarry Smith ierr = PetscMalloc4(am,PetscInt*,&uj_ptr,am,PetscInt*,&uj_lvl_ptr,am,PetscInt,&il,am,PetscInt,&jl);CHKERRQ(ierr); 103478910aadSHong Zhang for (i=0; i<am; i++){ 1035e60cf9a0SBarry Smith jl[i] = am; il[i] = 0; 103678910aadSHong Zhang } 103778910aadSHong Zhang 103878910aadSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 103978910aadSHong Zhang nlnk = am + 1; 104078910aadSHong Zhang ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 104178910aadSHong Zhang 104278910aadSHong Zhang /* initial FreeSpace size is fill*(ai[am]+1) */ 1043a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 104478910aadSHong Zhang current_space = free_space; 1045a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 104678910aadSHong Zhang current_space_lvl = free_space_lvl; 104778910aadSHong Zhang 104878910aadSHong Zhang for (k=0; k<am; k++){ /* for each active row k */ 104978910aadSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 105078910aadSHong Zhang nzk = 0; 105178910aadSHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 105278910aadSHong Zhang ncols_upper = 0; 105378910aadSHong Zhang cols = cols_lvl + am; 105478910aadSHong Zhang for (j=0; j<ncols; j++){ 105578910aadSHong Zhang i = rip[*(aj + ai[rip[k]] + j)]; 105678910aadSHong Zhang if (i >= k){ /* only take upper triangular entry */ 105778910aadSHong Zhang cols[ncols_upper] = i; 105878910aadSHong Zhang cols_lvl[ncols_upper] = -1; /* initialize level for nonzero entries */ 105978910aadSHong Zhang ncols_upper++; 106078910aadSHong Zhang } 106178910aadSHong Zhang } 106278910aadSHong Zhang ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 106378910aadSHong Zhang nzk += nlnk; 106478910aadSHong Zhang 106578910aadSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 106678910aadSHong Zhang prow = jl[k]; /* 1st pivot row */ 106778910aadSHong Zhang 106878910aadSHong Zhang while (prow < k){ 106978910aadSHong Zhang nextprow = jl[prow]; 107078910aadSHong Zhang 107178910aadSHong Zhang /* merge prow into k-th row */ 107278910aadSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 107378910aadSHong Zhang jmax = ui[prow+1]; 107478910aadSHong Zhang ncols = jmax-jmin; 107578910aadSHong Zhang i = jmin - ui[prow]; 107678910aadSHong Zhang cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 107778910aadSHong Zhang for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j); 10785a8e39fbSHong Zhang ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 107978910aadSHong Zhang nzk += nlnk; 108078910aadSHong Zhang 108178910aadSHong Zhang /* update il and jl for prow */ 108278910aadSHong Zhang if (jmin < jmax){ 108378910aadSHong Zhang il[prow] = jmin; 108478910aadSHong Zhang j = *cols; jl[prow] = jl[j]; jl[j] = prow; 108578910aadSHong Zhang } 108678910aadSHong Zhang prow = nextprow; 108778910aadSHong Zhang } 108878910aadSHong Zhang 108978910aadSHong Zhang /* if free space is not available, make more free space */ 109078910aadSHong Zhang if (current_space->local_remaining<nzk) { 109178910aadSHong Zhang i = am - k + 1; /* num of unfactored rows */ 109278910aadSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 1093a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 1094a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space_lvl);CHKERRQ(ierr); 109578910aadSHong Zhang reallocs++; 109678910aadSHong Zhang } 109778910aadSHong Zhang 109878910aadSHong Zhang /* copy data into free_space and free_space_lvl, then initialize lnk */ 109978910aadSHong Zhang ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 110078910aadSHong Zhang 110178910aadSHong Zhang /* add the k-th row into il and jl */ 110278910aadSHong Zhang if (nzk-1 > 0){ 110378910aadSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 110478910aadSHong Zhang jl[k] = jl[i]; jl[i] = k; 110578910aadSHong Zhang il[k] = ui[k] + 1; 110678910aadSHong Zhang } 110778910aadSHong Zhang uj_ptr[k] = current_space->array; 110878910aadSHong Zhang uj_lvl_ptr[k] = current_space_lvl->array; 110978910aadSHong Zhang 111078910aadSHong Zhang current_space->array += nzk; 111178910aadSHong Zhang current_space->local_used += nzk; 111278910aadSHong Zhang current_space->local_remaining -= nzk; 111378910aadSHong Zhang 111478910aadSHong Zhang current_space_lvl->array += nzk; 111578910aadSHong Zhang current_space_lvl->local_used += nzk; 111678910aadSHong Zhang current_space_lvl->local_remaining -= nzk; 111778910aadSHong Zhang 111878910aadSHong Zhang ui[k+1] = ui[k] + nzk; 111978910aadSHong Zhang } 112078910aadSHong Zhang 11216cf91177SBarry Smith #if defined(PETSC_USE_INFO) 112278910aadSHong Zhang if (ai[am] != 0) { 112378910aadSHong Zhang PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]); 1124ae15b995SBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 1125ae15b995SBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 1126ae15b995SBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 112778910aadSHong Zhang } else { 1128ae15b995SBarry Smith ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 112978910aadSHong Zhang } 113063ba0a88SBarry Smith #endif 113178910aadSHong Zhang 113278910aadSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 1133fca92195SBarry Smith ierr = PetscFree4(uj_ptr,uj_lvl_ptr,il,jl);CHKERRQ(ierr); 113478910aadSHong Zhang ierr = PetscFree(cols_lvl);CHKERRQ(ierr); 113578910aadSHong Zhang 113678910aadSHong Zhang /* destroy list of free space and other temporary array(s) */ 113778910aadSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 1138a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 113978910aadSHong Zhang ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1140a1a86e44SBarry Smith ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 114178910aadSHong Zhang 114278910aadSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 1143719d5645SBarry Smith B = fact; 1144ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 114578910aadSHong Zhang 114678910aadSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 114778910aadSHong Zhang b->singlemalloc = PETSC_FALSE; 1148e6b907acSBarry Smith b->free_a = PETSC_TRUE; 1149e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 115078910aadSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 115178910aadSHong Zhang b->j = uj; 115278910aadSHong Zhang b->i = ui; 115378910aadSHong Zhang b->diag = 0; 115478910aadSHong Zhang b->ilen = 0; 115578910aadSHong Zhang b->imax = 0; 115678910aadSHong Zhang b->row = perm; 115778910aadSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 115878910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 115978910aadSHong Zhang b->icol = perm; 116078910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 116178910aadSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 116252e6d16bSBarry Smith ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 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.) { 116878910aadSHong Zhang B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 116978910aadSHong Zhang } else { 117078910aadSHong Zhang B->info.fill_ratio_needed = 0.0; 117178910aadSHong Zhang } 117278910aadSHong Zhang if (perm_identity){ 11730a3c351aSHong Zhang B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 11740a3c351aSHong Zhang B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 117578910aadSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 117678910aadSHong Zhang } else { 1177719d5645SBarry Smith (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 117878910aadSHong Zhang } 1179c05c3958SHong Zhang PetscFunctionReturn(0); 1180c05c3958SHong Zhang } 1181c05c3958SHong Zhang 1182c05c3958SHong Zhang #undef __FUNCT__ 1183c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ" 11840481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 1185c05c3958SHong Zhang { 118678910aadSHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 118778910aadSHong Zhang Mat_SeqSBAIJ *b; 118878910aadSHong Zhang Mat B; 118978910aadSHong Zhang PetscErrorCode ierr; 119078910aadSHong Zhang PetscTruth perm_identity; 119178910aadSHong Zhang PetscReal fill = info->fill; 11925d0c19d7SBarry Smith const PetscInt *rip; 11935d0c19d7SBarry Smith PetscInt i,mbs=a->mbs,bs=A->rmap->bs,*ai=a->i,*aj=a->j,reallocs=0,prow; 119478910aadSHong Zhang PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 119578910aadSHong Zhang PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 1196a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 119778910aadSHong Zhang PetscBT lnkbt; 119878910aadSHong Zhang 1199c05c3958SHong Zhang PetscFunctionBegin; 12006ad2eaddSHong Zhang if (bs > 1) { /* convert to seqsbaij */ 12016ad2eaddSHong Zhang if (!a->sbaijMat){ 1202ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 12036ad2eaddSHong Zhang } 1204719d5645SBarry Smith (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */ 1205719d5645SBarry Smith ierr = MatCholeskyFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr); 12066ad2eaddSHong Zhang PetscFunctionReturn(0); 12076ad2eaddSHong Zhang } 12086ad2eaddSHong Zhang 120978910aadSHong Zhang /* check whether perm is the identity mapping */ 121078910aadSHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 1211e32f2f54SBarry Smith if (!perm_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported"); 121278910aadSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 121378910aadSHong Zhang 121478910aadSHong Zhang /* initialization */ 121578910aadSHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 1216e60cf9a0SBarry Smith ui[0] = 0; 121778910aadSHong Zhang 121878910aadSHong Zhang /* jl: linked list for storing indices of the pivot rows 121978910aadSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */ 1220fca92195SBarry Smith ierr = PetscMalloc4(mbs,PetscInt*,&ui_ptr,mbs,PetscInt,&il,mbs,PetscInt,&jl,mbs,PetscInt,&cols);CHKERRQ(ierr); 122178910aadSHong Zhang for (i=0; i<mbs; i++){ 1222e60cf9a0SBarry Smith jl[i] = mbs; il[i] = 0; 122378910aadSHong Zhang } 122478910aadSHong Zhang 122578910aadSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 122678910aadSHong Zhang nlnk = mbs + 1; 122778910aadSHong Zhang ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 122878910aadSHong Zhang 122978910aadSHong Zhang /* initial FreeSpace size is fill*(ai[mbs]+1) */ 1230a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr); 123178910aadSHong Zhang current_space = free_space; 123278910aadSHong Zhang 123378910aadSHong Zhang for (k=0; k<mbs; k++){ /* for each active row k */ 123478910aadSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 123578910aadSHong Zhang nzk = 0; 123678910aadSHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 1237e60cf9a0SBarry Smith ncols_upper = 0; 123878910aadSHong Zhang for (j=0; j<ncols; j++){ 123978910aadSHong Zhang i = rip[*(aj + ai[rip[k]] + j)]; 124078910aadSHong Zhang if (i >= k){ /* only take upper triangular entry */ 124178910aadSHong Zhang cols[ncols_upper] = i; 124278910aadSHong Zhang ncols_upper++; 124378910aadSHong Zhang } 124478910aadSHong Zhang } 124578910aadSHong Zhang ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 124678910aadSHong Zhang nzk += nlnk; 124778910aadSHong Zhang 124878910aadSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 124978910aadSHong Zhang prow = jl[k]; /* 1st pivot row */ 125078910aadSHong Zhang 125178910aadSHong Zhang while (prow < k){ 125278910aadSHong Zhang nextprow = jl[prow]; 125378910aadSHong Zhang /* merge prow into k-th row */ 125478910aadSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */ 125578910aadSHong Zhang jmax = ui[prow+1]; 125678910aadSHong Zhang ncols = jmax-jmin; 125778910aadSHong Zhang uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */ 12585a8e39fbSHong Zhang ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 125978910aadSHong Zhang nzk += nlnk; 126078910aadSHong Zhang 126178910aadSHong Zhang /* update il and jl for prow */ 126278910aadSHong Zhang if (jmin < jmax){ 126378910aadSHong Zhang il[prow] = jmin; 126478910aadSHong Zhang j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 126578910aadSHong Zhang } 126678910aadSHong Zhang prow = nextprow; 126778910aadSHong Zhang } 126878910aadSHong Zhang 126978910aadSHong Zhang /* if free space is not available, make more free space */ 127078910aadSHong Zhang if (current_space->local_remaining<nzk) { 127178910aadSHong Zhang i = mbs - k + 1; /* num of unfactored rows */ 127278910aadSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 1273a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 127478910aadSHong Zhang reallocs++; 127578910aadSHong Zhang } 127678910aadSHong Zhang 127778910aadSHong Zhang /* copy data into free space, then initialize lnk */ 127878910aadSHong Zhang ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 127978910aadSHong Zhang 128078910aadSHong Zhang /* add the k-th row into il and jl */ 128178910aadSHong Zhang if (nzk-1 > 0){ 128278910aadSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ 128378910aadSHong Zhang jl[k] = jl[i]; jl[i] = k; 128478910aadSHong Zhang il[k] = ui[k] + 1; 128578910aadSHong Zhang } 128678910aadSHong Zhang ui_ptr[k] = current_space->array; 128778910aadSHong Zhang current_space->array += nzk; 128878910aadSHong Zhang current_space->local_used += nzk; 128978910aadSHong Zhang current_space->local_remaining -= nzk; 129078910aadSHong Zhang 129178910aadSHong Zhang ui[k+1] = ui[k] + nzk; 129278910aadSHong Zhang } 129378910aadSHong Zhang 12946cf91177SBarry Smith #if defined(PETSC_USE_INFO) 129575567043SBarry Smith if (ai[mbs] != 0.) { 129678910aadSHong Zhang PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 1297ae15b995SBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 1298ae15b995SBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 1299ae15b995SBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 130078910aadSHong Zhang } else { 1301ae15b995SBarry Smith ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 130278910aadSHong Zhang } 130363ba0a88SBarry Smith #endif 130478910aadSHong Zhang 130578910aadSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 1306fca92195SBarry Smith ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr); 130778910aadSHong Zhang 130878910aadSHong Zhang /* destroy list of free space and other temporary array(s) */ 130978910aadSHong Zhang ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 1310a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 131178910aadSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 131278910aadSHong Zhang 131378910aadSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 1314719d5645SBarry Smith B = fact; 1315ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 131678910aadSHong Zhang 131778910aadSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 131878910aadSHong Zhang b->singlemalloc = PETSC_FALSE; 1319e6b907acSBarry Smith b->free_a = PETSC_TRUE; 1320e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 132178910aadSHong Zhang ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 132278910aadSHong Zhang b->j = uj; 132378910aadSHong Zhang b->i = ui; 132478910aadSHong Zhang b->diag = 0; 132578910aadSHong Zhang b->ilen = 0; 132678910aadSHong Zhang b->imax = 0; 132778910aadSHong Zhang b->row = perm; 132878910aadSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 132978910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 133078910aadSHong Zhang b->icol = perm; 133178910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 133278910aadSHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 133352e6d16bSBarry Smith ierr = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 133478910aadSHong Zhang b->maxnz = b->nz = ui[mbs]; 133578910aadSHong Zhang 133678910aadSHong Zhang B->info.factor_mallocs = reallocs; 133778910aadSHong Zhang B->info.fill_ratio_given = fill; 133875567043SBarry Smith if (ai[mbs] != 0.) { 133978910aadSHong Zhang B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 134078910aadSHong Zhang } else { 134178910aadSHong Zhang B->info.fill_ratio_needed = 0.0; 134278910aadSHong Zhang } 134378910aadSHong Zhang if (perm_identity){ 13446ad2eaddSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 134578910aadSHong Zhang } else { 13466ad2eaddSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 134778910aadSHong Zhang } 1348c05c3958SHong Zhang PetscFunctionReturn(0); 1349c05c3958SHong Zhang } 1350c8342467SHong Zhang 13512efa7f71SHong Zhang #undef __FUNCT__ 13524dd39f65SShri Abhyankar #define __FUNCT__ "MatSolve_SeqBAIJ_N_NaturalOrdering" 13534dd39f65SShri Abhyankar PetscErrorCode MatSolve_SeqBAIJ_N_NaturalOrdering(Mat A,Vec bb,Vec xx) 13541a83e813SShri Abhyankar { 13551a83e813SShri Abhyankar Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 13561a83e813SShri Abhyankar PetscErrorCode ierr; 13571a83e813SShri Abhyankar const PetscInt *ai=a->i,*aj=a->j,*adiag=a->diag,*vi; 13581a83e813SShri Abhyankar PetscInt i,k,n=a->mbs; 13591a83e813SShri Abhyankar PetscInt nz,bs=A->rmap->bs,bs2=a->bs2; 13601a83e813SShri Abhyankar MatScalar *aa=a->a,*v; 13611a83e813SShri Abhyankar PetscScalar *x,*b,*s,*t,*ls; 13621a83e813SShri Abhyankar 13631a83e813SShri Abhyankar PetscFunctionBegin; 13641a83e813SShri Abhyankar ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 13651a83e813SShri Abhyankar ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 13661a83e813SShri Abhyankar t = a->solve_work; 13671a83e813SShri Abhyankar 13681a83e813SShri Abhyankar /* forward solve the lower triangular */ 13691a83e813SShri Abhyankar ierr = PetscMemcpy(t,b,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy 1st block of b to t */ 13701a83e813SShri Abhyankar 13711a83e813SShri Abhyankar for (i=1; i<n; i++) { 13721a83e813SShri Abhyankar v = aa + bs2*ai[i]; 13731a83e813SShri Abhyankar vi = aj + ai[i]; 13741a83e813SShri Abhyankar nz = ai[i+1] - ai[i]; 13751a83e813SShri Abhyankar s = t + bs*i; 13761a83e813SShri Abhyankar ierr = PetscMemcpy(s,b+bs*i,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy i_th block of b to t */ 13771a83e813SShri Abhyankar for(k=0;k<nz;k++){ 13781a83e813SShri Abhyankar Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[k]); 13791a83e813SShri Abhyankar v += bs2; 13801a83e813SShri Abhyankar } 13811a83e813SShri Abhyankar } 13821a83e813SShri Abhyankar 13831a83e813SShri Abhyankar /* backward solve the upper triangular */ 13841a83e813SShri Abhyankar ls = a->solve_work + A->cmap->n; 13851a83e813SShri Abhyankar for (i=n-1; i>=0; i--){ 13861a83e813SShri Abhyankar v = aa + bs2*(adiag[i+1]+1); 13871a83e813SShri Abhyankar vi = aj + adiag[i+1]+1; 13881a83e813SShri Abhyankar nz = adiag[i] - adiag[i+1]-1; 13891a83e813SShri Abhyankar ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 13901a83e813SShri Abhyankar for(k=0;k<nz;k++){ 13911a83e813SShri Abhyankar Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[k]); 13921a83e813SShri Abhyankar v += bs2; 13931a83e813SShri Abhyankar } 13941a83e813SShri Abhyankar Kernel_w_gets_A_times_v(bs,ls,aa+bs2*adiag[i],t+i*bs); /* *inv(diagonal[i]) */ 13951a83e813SShri Abhyankar ierr = PetscMemcpy(x+i*bs,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 13961a83e813SShri Abhyankar } 13971a83e813SShri Abhyankar 13981a83e813SShri Abhyankar ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 13991a83e813SShri Abhyankar ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 14001a83e813SShri Abhyankar ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr); 14011a83e813SShri Abhyankar PetscFunctionReturn(0); 14021a83e813SShri Abhyankar } 14031a83e813SShri Abhyankar 14042efa7f71SHong Zhang #undef __FUNCT__ 14054dd39f65SShri Abhyankar #define __FUNCT__ "MatSolve_SeqBAIJ_N" 14064dd39f65SShri Abhyankar PetscErrorCode MatSolve_SeqBAIJ_N(Mat A,Vec bb,Vec xx) 140735aa4fcfSShri Abhyankar { 140835aa4fcfSShri Abhyankar Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 140935aa4fcfSShri Abhyankar IS iscol=a->col,isrow=a->row; 141035aa4fcfSShri Abhyankar PetscErrorCode ierr; 141135aa4fcfSShri Abhyankar const PetscInt *r,*c,*rout,*cout,*ai=a->i,*aj=a->j,*adiag=a->diag,*vi; 141235aa4fcfSShri Abhyankar PetscInt i,m,n=a->mbs; 141335aa4fcfSShri Abhyankar PetscInt nz,bs=A->rmap->bs,bs2=a->bs2; 141435aa4fcfSShri Abhyankar MatScalar *aa=a->a,*v; 141535aa4fcfSShri Abhyankar PetscScalar *x,*b,*s,*t,*ls; 141635aa4fcfSShri Abhyankar 141735aa4fcfSShri Abhyankar PetscFunctionBegin; 141835aa4fcfSShri Abhyankar ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 141935aa4fcfSShri Abhyankar ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 142035aa4fcfSShri Abhyankar t = a->solve_work; 142135aa4fcfSShri Abhyankar 142235aa4fcfSShri Abhyankar ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 142335aa4fcfSShri Abhyankar ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 142435aa4fcfSShri Abhyankar 142535aa4fcfSShri Abhyankar /* forward solve the lower triangular */ 142635aa4fcfSShri Abhyankar ierr = PetscMemcpy(t,b+bs*r[0],bs*sizeof(PetscScalar));CHKERRQ(ierr); 142735aa4fcfSShri Abhyankar for (i=1; i<n; i++) { 142835aa4fcfSShri Abhyankar v = aa + bs2*ai[i]; 142935aa4fcfSShri Abhyankar vi = aj + ai[i]; 143035aa4fcfSShri Abhyankar nz = ai[i+1] - ai[i]; 143135aa4fcfSShri Abhyankar s = t + bs*i; 143235aa4fcfSShri Abhyankar ierr = PetscMemcpy(s,b+bs*r[i],bs*sizeof(PetscScalar));CHKERRQ(ierr); 143335aa4fcfSShri Abhyankar for(m=0;m<nz;m++){ 143435aa4fcfSShri Abhyankar Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[m]); 143535aa4fcfSShri Abhyankar v += bs2; 143635aa4fcfSShri Abhyankar } 143735aa4fcfSShri Abhyankar } 143835aa4fcfSShri Abhyankar 143935aa4fcfSShri Abhyankar /* backward solve the upper triangular */ 144035aa4fcfSShri Abhyankar ls = a->solve_work + A->cmap->n; 144135aa4fcfSShri Abhyankar for (i=n-1; i>=0; i--){ 144235aa4fcfSShri Abhyankar v = aa + bs2*(adiag[i+1]+1); 144335aa4fcfSShri Abhyankar vi = aj + adiag[i+1]+1; 144435aa4fcfSShri Abhyankar nz = adiag[i] - adiag[i+1] - 1; 144535aa4fcfSShri Abhyankar ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 144635aa4fcfSShri Abhyankar for(m=0;m<nz;m++){ 144735aa4fcfSShri Abhyankar Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[m]); 144835aa4fcfSShri Abhyankar v += bs2; 144935aa4fcfSShri Abhyankar } 145035aa4fcfSShri Abhyankar Kernel_w_gets_A_times_v(bs,ls,v,t+i*bs); /* *inv(diagonal[i]) */ 145135aa4fcfSShri Abhyankar ierr = PetscMemcpy(x + bs*c[i],t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 145235aa4fcfSShri Abhyankar } 145335aa4fcfSShri Abhyankar ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 145435aa4fcfSShri Abhyankar ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 145535aa4fcfSShri Abhyankar ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 145635aa4fcfSShri Abhyankar ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 145735aa4fcfSShri Abhyankar ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr); 145835aa4fcfSShri Abhyankar PetscFunctionReturn(0); 145935aa4fcfSShri Abhyankar } 146035aa4fcfSShri Abhyankar 146135aa4fcfSShri Abhyankar #undef __FUNCT__ 14622efa7f71SHong Zhang #define __FUNCT__ "BlockAbs_privat" 146391d91c02SMatthew Knepley PetscErrorCode BlockAbs_private(PetscInt nbs,PetscInt bs2,PetscScalar *blockarray,PetscReal *absarray) 14642efa7f71SHong Zhang { 14652efa7f71SHong Zhang PetscErrorCode ierr; 14662efa7f71SHong Zhang PetscInt i,j; 14672efa7f71SHong Zhang PetscFunctionBegin; 14682efa7f71SHong Zhang ierr = PetscMemzero(absarray,(nbs+1)*sizeof(PetscReal));CHKERRQ(ierr); 14692efa7f71SHong Zhang for (i=0; i<nbs; i++){ 14702efa7f71SHong Zhang for (j=0; j<bs2; j++){ 14712efa7f71SHong Zhang if (absarray[i] < PetscAbsScalar(blockarray[i*nbs+j])) absarray[i] = PetscAbsScalar(blockarray[i*nbs+j]); 14722efa7f71SHong Zhang } 14732efa7f71SHong Zhang } 14742efa7f71SHong Zhang PetscFunctionReturn(0); 14752efa7f71SHong Zhang } 14762efa7f71SHong Zhang 1477c8342467SHong Zhang #undef __FUNCT__ 1478c8342467SHong Zhang #define __FUNCT__ "MatILUDTFactor_SeqBAIJ" 1479fe97e370SBarry Smith /* 1480fe97e370SBarry Smith This needs to be renamed and called by the regular MatILUFactor_SeqBAIJ when drop tolerance is used 1481fe97e370SBarry Smith */ 1482c8342467SHong Zhang PetscErrorCode MatILUDTFactor_SeqBAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact) 1483c8342467SHong Zhang { 14842efa7f71SHong Zhang Mat B = *fact; 14852efa7f71SHong Zhang Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b; 14862efa7f71SHong Zhang IS isicol; 14872efa7f71SHong Zhang PetscErrorCode ierr; 14882efa7f71SHong Zhang const PetscInt *r,*ic; 14892efa7f71SHong Zhang PetscInt i,mbs=a->mbs,bs=A->rmap->bs,bs2=a->bs2,*ai=a->i,*aj=a->j,*ajtmp,*adiag; 14902efa7f71SHong Zhang PetscInt *bi,*bj,*bdiag; 14912efa7f71SHong Zhang 14922efa7f71SHong Zhang PetscInt row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au; 14932efa7f71SHong Zhang PetscInt nlnk,*lnk; 14942efa7f71SHong Zhang PetscBT lnkbt; 14952efa7f71SHong Zhang PetscTruth row_identity,icol_identity,both_identity; 14962efa7f71SHong Zhang MatScalar *aatmp,*pv,*batmp,*ba,*rtmp,*pc,*multiplier,*vtmp; 14972efa7f71SHong Zhang const PetscInt *ics; 14982efa7f71SHong Zhang PetscInt j,nz,*pj,*bjtmp,k,ncut,*jtmp; 14992efa7f71SHong Zhang 1500182b8fbaSHong Zhang PetscReal dt=info->dt; /* shift=info->shiftamount; */ 15012efa7f71SHong Zhang PetscInt nnz_max; 15022efa7f71SHong Zhang PetscTruth missing; 1503021822bcSHong Zhang PetscReal *vtmp_abs; 150497ef94ebSSatish Balay MatScalar *v_work; 150597ef94ebSSatish Balay PetscInt *v_pivots; 15062efa7f71SHong Zhang 1507c8342467SHong Zhang PetscFunctionBegin; 15082efa7f71SHong Zhang /* ------- symbolic factorization, can be reused ---------*/ 15092efa7f71SHong Zhang ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr); 1510e32f2f54SBarry Smith if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i); 15112efa7f71SHong Zhang adiag=a->diag; 15122efa7f71SHong Zhang 15132efa7f71SHong Zhang ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 15142efa7f71SHong Zhang 15152efa7f71SHong Zhang /* bdiag is location of diagonal in factor */ 15162efa7f71SHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 15172efa7f71SHong Zhang 15182efa7f71SHong Zhang /* allocate row pointers bi */ 15196bce7ff8SHong Zhang ierr = PetscMalloc((2*mbs+2)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 15202efa7f71SHong Zhang 15212efa7f71SHong Zhang /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */ 15222efa7f71SHong Zhang dtcount = (PetscInt)info->dtcount; 15236bce7ff8SHong Zhang if (dtcount > mbs-1) dtcount = mbs-1; 15246bce7ff8SHong Zhang nnz_max = ai[mbs]+2*mbs*dtcount +2; 15256da515a1SHong Zhang /* printf("MatILUDTFactor_SeqBAIJ, bs %d, ai[mbs] %d, nnz_max %d, dtcount %d\n",bs,ai[mbs],nnz_max,dtcount); */ 15262efa7f71SHong Zhang ierr = PetscMalloc(nnz_max*sizeof(PetscInt),&bj);CHKERRQ(ierr); 15276bce7ff8SHong Zhang nnz_max = nnz_max*bs2; 15282efa7f71SHong Zhang ierr = PetscMalloc(nnz_max*sizeof(MatScalar),&ba);CHKERRQ(ierr); 15292efa7f71SHong Zhang 15302efa7f71SHong Zhang /* put together the new matrix */ 15312efa7f71SHong Zhang ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 15322efa7f71SHong Zhang ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr); 15332efa7f71SHong Zhang b = (Mat_SeqBAIJ*)(B)->data; 15342efa7f71SHong Zhang b->free_a = PETSC_TRUE; 15352efa7f71SHong Zhang b->free_ij = PETSC_TRUE; 15362efa7f71SHong Zhang b->singlemalloc = PETSC_FALSE; 15372efa7f71SHong Zhang b->a = ba; 15382efa7f71SHong Zhang b->j = bj; 15392efa7f71SHong Zhang b->i = bi; 15402efa7f71SHong Zhang b->diag = bdiag; 15412efa7f71SHong Zhang b->ilen = 0; 15422efa7f71SHong Zhang b->imax = 0; 15432efa7f71SHong Zhang b->row = isrow; 15442efa7f71SHong Zhang b->col = iscol; 15452efa7f71SHong Zhang ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 15462efa7f71SHong Zhang ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 15472efa7f71SHong Zhang b->icol = isicol; 15482efa7f71SHong Zhang ierr = PetscMalloc((bs*(mbs+1))*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 15492efa7f71SHong Zhang 15502efa7f71SHong Zhang ierr = PetscLogObjectMemory(B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 15512efa7f71SHong Zhang b->maxnz = nnz_max/bs2; 15522efa7f71SHong Zhang 1553d5f3da31SBarry Smith (B)->factortype = MAT_FACTOR_ILUDT; 15542efa7f71SHong Zhang (B)->info.factor_mallocs = 0; 15552efa7f71SHong Zhang (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)(ai[mbs]*bs2)); 15562efa7f71SHong Zhang CHKMEMQ; 15572efa7f71SHong Zhang /* ------- end of symbolic factorization ---------*/ 15582efa7f71SHong Zhang ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 15592efa7f71SHong Zhang ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 15602efa7f71SHong Zhang ics = ic; 15612efa7f71SHong Zhang 15622efa7f71SHong Zhang /* linked list for storing column indices of the active row */ 15632efa7f71SHong Zhang nlnk = mbs + 1; 15642efa7f71SHong Zhang ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 15652efa7f71SHong Zhang 15662efa7f71SHong Zhang /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */ 1567fca92195SBarry Smith ierr = PetscMalloc2(mbs,PetscInt,&im,mbs,PetscInt,&jtmp);CHKERRQ(ierr); 15682efa7f71SHong Zhang /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */ 1569fca92195SBarry Smith ierr = PetscMalloc2(mbs*bs2,MatScalar,&rtmp,mbs*bs2,MatScalar,&vtmp);CHKERRQ(ierr); 1570021822bcSHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(PetscReal),&vtmp_abs);CHKERRQ(ierr); 1571fca92195SBarry Smith ierr = PetscMalloc3(bs,MatScalar,&v_work,bs2,MatScalar,&multiplier,bs,PetscInt,&v_pivots);CHKERRQ(ierr); 15722efa7f71SHong Zhang 1573e60cf9a0SBarry Smith bi[0] = 0; 15742efa7f71SHong Zhang bdiag[0] = (nnz_max/bs2)-1; /* location of diagonal in factor B */ 15756bce7ff8SHong Zhang bi[2*mbs+1] = bdiag[0]+1; /* endof bj and ba array */ 15762efa7f71SHong Zhang for (i=0; i<mbs; i++) { 15772efa7f71SHong Zhang /* copy initial fill into linked list */ 15782efa7f71SHong Zhang nzi = 0; /* nonzeros for active row i */ 15792efa7f71SHong Zhang nzi = ai[r[i]+1] - ai[r[i]]; 1580e32f2f54SBarry Smith if (!nzi) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 15812efa7f71SHong Zhang nzi_al = adiag[r[i]] - ai[r[i]]; 15822efa7f71SHong Zhang nzi_au = ai[r[i]+1] - adiag[r[i]] -1; 15836da515a1SHong Zhang /* printf("row %d, nzi_al/au %d %d\n",i,nzi_al,nzi_au); */ 15842efa7f71SHong Zhang 15852efa7f71SHong Zhang /* load in initial unfactored row */ 15862efa7f71SHong Zhang ajtmp = aj + ai[r[i]]; 15872efa7f71SHong Zhang ierr = PetscLLAddPerm(nzi,ajtmp,ic,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1588fca92195SBarry Smith ierr = PetscMemzero(rtmp,mbs*bs2*sizeof(PetscScalar));CHKERRQ(ierr); 15892efa7f71SHong Zhang aatmp = a->a + bs2*ai[r[i]]; 15902efa7f71SHong Zhang for (j=0; j<nzi; j++) { 15912efa7f71SHong Zhang ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],aatmp+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 15922efa7f71SHong Zhang } 15932efa7f71SHong Zhang 15942efa7f71SHong Zhang /* add pivot rows into linked list */ 15952efa7f71SHong Zhang row = lnk[mbs]; 15962efa7f71SHong Zhang while (row < i) { 15972efa7f71SHong Zhang nzi_bl = bi[row+1] - bi[row] + 1; 15982efa7f71SHong Zhang bjtmp = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */ 15992efa7f71SHong Zhang ierr = PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);CHKERRQ(ierr); 16002efa7f71SHong Zhang nzi += nlnk; 16012efa7f71SHong Zhang row = lnk[row]; 16022efa7f71SHong Zhang } 16032efa7f71SHong Zhang 16042efa7f71SHong Zhang /* copy data from lnk into jtmp, then initialize lnk */ 16052efa7f71SHong Zhang ierr = PetscLLClean(mbs,mbs,nzi,lnk,jtmp,lnkbt);CHKERRQ(ierr); 16062efa7f71SHong Zhang 16072efa7f71SHong Zhang /* numerical factorization */ 16082efa7f71SHong Zhang bjtmp = jtmp; 16092efa7f71SHong Zhang row = *bjtmp++; /* 1st pivot row */ 16102efa7f71SHong Zhang 16112efa7f71SHong Zhang while (row < i) { 16122efa7f71SHong Zhang pc = rtmp + bs2*row; 16132efa7f71SHong Zhang pv = ba + bs2*bdiag[row]; /* inv(diag) of the pivot row */ 16142efa7f71SHong Zhang Kernel_A_gets_A_times_B(bs,pc,pv,multiplier); /* pc= multiplier = pc*inv(diag[row]) */ 16152efa7f71SHong Zhang ierr = BlockAbs_private(1,bs2,pc,vtmp_abs);CHKERRQ(ierr); 16162efa7f71SHong Zhang if (vtmp_abs[0] > dt){ /* apply tolerance dropping rule */ 16172efa7f71SHong Zhang pj = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */ 16182efa7f71SHong Zhang pv = ba + bs2*(bdiag[row+1] + 1); 16192efa7f71SHong Zhang nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */ 16202efa7f71SHong Zhang for (j=0; j<nz; j++){ 16212efa7f71SHong Zhang Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); 16222efa7f71SHong Zhang } 16236da515a1SHong Zhang /* ierr = PetscLogFlops(bslog*(nz+1.0)-bs);CHKERRQ(ierr); */ 16242efa7f71SHong Zhang } 16252efa7f71SHong Zhang row = *bjtmp++; 16262efa7f71SHong Zhang } 16272efa7f71SHong Zhang 16282efa7f71SHong Zhang /* copy sparse rtmp into contiguous vtmp; separate L and U part */ 16292efa7f71SHong Zhang nzi_bl = 0; j = 0; 16302efa7f71SHong Zhang while (jtmp[j] < i){ /* L-part. Note: jtmp is sorted */ 16312efa7f71SHong Zhang ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 16322efa7f71SHong Zhang nzi_bl++; j++; 16332efa7f71SHong Zhang } 16342efa7f71SHong Zhang nzi_bu = nzi - nzi_bl -1; 16356da515a1SHong Zhang /* printf("nzi %d, nzi_bl %d, nzi_bu %d\n",nzi,nzi_bl,nzi_bu); */ 16362efa7f71SHong Zhang 16372efa7f71SHong Zhang while (j < nzi){ /* U-part */ 16382efa7f71SHong Zhang ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 16392efa7f71SHong Zhang /* 16402efa7f71SHong Zhang printf(" col %d: ",jtmp[j]); 16412efa7f71SHong Zhang for (j1=0; j1<bs2; j1++) printf(" %g",*(vtmp+bs2*j+j1)); 16422efa7f71SHong Zhang printf(" \n"); 16432efa7f71SHong Zhang */ 16442efa7f71SHong Zhang j++; 16452efa7f71SHong Zhang } 16462efa7f71SHong Zhang 16472efa7f71SHong Zhang ierr = BlockAbs_private(nzi,bs2,vtmp,vtmp_abs);CHKERRQ(ierr); 16482efa7f71SHong Zhang /* 16492efa7f71SHong Zhang printf(" row %d, nzi %d, vtmp_abs\n",i,nzi); 16502efa7f71SHong Zhang for (j1=0; j1<nzi; j1++) printf(" (%d %g),",jtmp[j1],vtmp_abs[j1]); 16512efa7f71SHong Zhang printf(" \n"); 16522efa7f71SHong 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){ 1658021822bcSHong Zhang ierr = PetscSortSplitReal(ncut,nzi_bl,vtmp_abs,jtmp);CHKERRQ(ierr); 16592efa7f71SHong Zhang ierr = PetscSortIntWithScalarArray(ncut,jtmp,vtmp);CHKERRQ(ierr); 16602efa7f71SHong Zhang } else { 16612efa7f71SHong Zhang ncut = nzi_bl; 16622efa7f71SHong Zhang } 16632efa7f71SHong Zhang for (j=0; j<ncut; j++){ 16642efa7f71SHong Zhang bjtmp[j] = jtmp[j]; 16652efa7f71SHong Zhang ierr = PetscMemcpy(batmp+bs2*j,rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 16662efa7f71SHong Zhang /* 16672efa7f71SHong Zhang printf(" col %d: ",bjtmp[j]); 16682efa7f71SHong Zhang for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*j+j1)); 16692efa7f71SHong Zhang printf("\n"); 16702efa7f71SHong Zhang */ 16712efa7f71SHong Zhang } 16722efa7f71SHong Zhang bi[i+1] = bi[i] + ncut; 16732efa7f71SHong Zhang nzi = ncut + 1; 16742efa7f71SHong Zhang 16752efa7f71SHong Zhang /* apply level dropping rule to U part */ 16762efa7f71SHong Zhang ncut = nzi_au + dtcount; 16772efa7f71SHong Zhang if (ncut < nzi_bu){ 1678021822bcSHong Zhang ierr = PetscSortSplitReal(ncut,nzi_bu,vtmp_abs+nzi_bl+1,jtmp+nzi_bl+1);CHKERRQ(ierr); 16792efa7f71SHong Zhang ierr = PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);CHKERRQ(ierr); 16802efa7f71SHong Zhang } else { 16812efa7f71SHong Zhang ncut = nzi_bu; 16822efa7f71SHong Zhang } 16832efa7f71SHong Zhang nzi += ncut; 16842efa7f71SHong Zhang 16852efa7f71SHong Zhang /* mark bdiagonal */ 16862efa7f71SHong Zhang bdiag[i+1] = bdiag[i] - (ncut + 1); 16876bce7ff8SHong Zhang bi[2*mbs - i] = bi[2*mbs - i +1] - (ncut + 1); 16886bce7ff8SHong Zhang 16892efa7f71SHong Zhang bjtmp = bj + bdiag[i]; 16902efa7f71SHong Zhang batmp = ba + bs2*bdiag[i]; 16912efa7f71SHong Zhang ierr = PetscMemcpy(batmp,rtmp+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 16922efa7f71SHong Zhang *bjtmp = i; 16932efa7f71SHong Zhang /* 16942efa7f71SHong Zhang printf(" diag %d: ",*bjtmp); 16952efa7f71SHong Zhang for (j=0; j<bs2; j++){ 16962efa7f71SHong Zhang printf(" %g,",batmp[j]); 16972efa7f71SHong Zhang } 16982efa7f71SHong Zhang printf("\n"); 16992efa7f71SHong Zhang */ 17002efa7f71SHong Zhang bjtmp = bj + bdiag[i+1]+1; 17012efa7f71SHong Zhang batmp = ba + (bdiag[i+1]+1)*bs2; 17022efa7f71SHong Zhang 17032efa7f71SHong Zhang for (k=0; k<ncut; k++){ 17042efa7f71SHong Zhang bjtmp[k] = jtmp[nzi_bl+1+k]; 17052efa7f71SHong Zhang ierr = PetscMemcpy(batmp+bs2*k,rtmp+bs2*bjtmp[k],bs2*sizeof(MatScalar));CHKERRQ(ierr); 17062efa7f71SHong Zhang /* 17072efa7f71SHong Zhang printf(" col %d:",bjtmp[k]); 17082efa7f71SHong Zhang for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*k+j1)); 17092efa7f71SHong Zhang printf("\n"); 17102efa7f71SHong Zhang */ 17112efa7f71SHong Zhang } 17122efa7f71SHong Zhang 17132efa7f71SHong Zhang im[i] = nzi; /* used by PetscLLAddSortedLU() */ 17142efa7f71SHong Zhang 17152efa7f71SHong Zhang /* invert diagonal block for simplier triangular solves - add shift??? */ 17162efa7f71SHong Zhang batmp = ba + bs2*bdiag[i]; 17172efa7f71SHong Zhang ierr = Kernel_A_gets_inverse_A(bs,batmp,v_pivots,v_work);CHKERRQ(ierr); 17182efa7f71SHong Zhang } /* for (i=0; i<mbs; i++) */ 1719fca92195SBarry Smith ierr = PetscFree3(v_work,multiplier,v_pivots);CHKERRQ(ierr); 17202efa7f71SHong Zhang 17216da515a1SHong Zhang /* printf("end of L %d, beginning of U %d\n",bi[mbs],bdiag[mbs]); */ 1722e32f2f54SBarry Smith if (bi[mbs] >= bdiag[mbs]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"end of L array %d cannot >= the beginning of U array %d",bi[mbs],bdiag[mbs]); 17232efa7f71SHong Zhang 17242efa7f71SHong Zhang ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 17252efa7f71SHong Zhang ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 17262efa7f71SHong Zhang 17272efa7f71SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 17282efa7f71SHong Zhang 1729fca92195SBarry Smith ierr = PetscFree2(im,jtmp);CHKERRQ(ierr); 1730fca92195SBarry Smith ierr = PetscFree2(rtmp,vtmp);CHKERRQ(ierr); 17312efa7f71SHong Zhang 17322efa7f71SHong Zhang ierr = PetscLogFlops(bs2*B->cmap->n);CHKERRQ(ierr); 17332efa7f71SHong Zhang b->maxnz = b->nz = bi[mbs] + bdiag[0] - bdiag[mbs]; 17342efa7f71SHong Zhang 17352efa7f71SHong Zhang ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 17362efa7f71SHong Zhang ierr = ISIdentity(isicol,&icol_identity);CHKERRQ(ierr); 17372efa7f71SHong Zhang both_identity = (PetscTruth) (row_identity && icol_identity); 17382efa7f71SHong Zhang if (row_identity && icol_identity) { 17394dd39f65SShri Abhyankar B->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering; 17402efa7f71SHong Zhang } else { 17414dd39f65SShri Abhyankar B->ops->solve = MatSolve_SeqBAIJ_N; 17422efa7f71SHong Zhang } 17432efa7f71SHong Zhang 17442efa7f71SHong Zhang B->ops->solveadd = 0; 17452efa7f71SHong Zhang B->ops->solvetranspose = 0; 17462efa7f71SHong Zhang B->ops->solvetransposeadd = 0; 17472efa7f71SHong Zhang B->ops->matsolve = 0; 17482efa7f71SHong Zhang B->assembled = PETSC_TRUE; 17492efa7f71SHong Zhang B->preallocated = PETSC_TRUE; 1750c8342467SHong Zhang PetscFunctionReturn(0); 1751c8342467SHong Zhang } 1752