1be1d678aSKris Buschelman 25d517e7eSBarry Smith /* 3ec3a8b7bSBarry Smith Factorization code for BAIJ format. 45d517e7eSBarry Smith */ 5c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h> 606873bf2SBarry Smith #include <petsc-private/kernels/blockinvert.h> 7ec3a8b7bSBarry Smith 8b588c5a2SHong Zhang #undef __FUNCT__ 94dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2" 104dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat B,Mat A,const MatFactorInfo *info) 11b588c5a2SHong Zhang { 12b588c5a2SHong Zhang Mat C =B; 13b588c5a2SHong Zhang Mat_SeqBAIJ *a =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data; 14b588c5a2SHong Zhang IS isrow = b->row,isicol = b->icol; 15b588c5a2SHong Zhang PetscErrorCode ierr; 165a586d82SBarry Smith const PetscInt *r,*ic; 17bbd65245SShri Abhyankar PetscInt i,j,k,nz,nzL,row,*pj; 18bbd65245SShri Abhyankar const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,bs2=a->bs2; 19bbd65245SShri Abhyankar const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag; 20bbd65245SShri Abhyankar MatScalar *rtmp,*pc,*mwork,*pv; 21bbd65245SShri Abhyankar MatScalar *aa=a->a,*v; 22bbd65245SShri Abhyankar PetscInt flg; 23182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 24b588c5a2SHong Zhang 25b588c5a2SHong Zhang PetscFunctionBegin; 26b588c5a2SHong Zhang ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 27b588c5a2SHong Zhang ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 28b588c5a2SHong Zhang 294c000e2eSHong Zhang /* generate work space needed by the factorization */ 30dcca6d9dSJed Brown ierr = PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);CHKERRQ(ierr); 314c000e2eSHong Zhang ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 32b588c5a2SHong Zhang 33b588c5a2SHong Zhang for (i=0; i<n; i++) { 34b588c5a2SHong Zhang /* zero rtmp */ 35b588c5a2SHong Zhang /* L part */ 36b588c5a2SHong Zhang nz = bi[i+1] - bi[i]; 37b588c5a2SHong Zhang bjtmp = bj + bi[i]; 38b588c5a2SHong Zhang for (j=0; j<nz; j++) { 39b588c5a2SHong Zhang ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 40b588c5a2SHong Zhang } 41b588c5a2SHong Zhang 42b588c5a2SHong Zhang /* U part */ 430c4413a7SShri Abhyankar nz = bdiag[i] - bdiag[i+1]; 440c4413a7SShri Abhyankar bjtmp = bj + bdiag[i+1]+1; 450c4413a7SShri Abhyankar for (j=0; j<nz; j++) { 460c4413a7SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 470c4413a7SShri Abhyankar } 480c4413a7SShri Abhyankar 490c4413a7SShri Abhyankar /* load in initial (unfactored row) */ 500c4413a7SShri Abhyankar nz = ai[r[i]+1] - ai[r[i]]; 510c4413a7SShri Abhyankar ajtmp = aj + ai[r[i]]; 520c4413a7SShri Abhyankar v = aa + bs2*ai[r[i]]; 530c4413a7SShri Abhyankar for (j=0; j<nz; j++) { 540c4413a7SShri Abhyankar ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 550c4413a7SShri Abhyankar } 560c4413a7SShri Abhyankar 570c4413a7SShri Abhyankar /* elimination */ 580c4413a7SShri Abhyankar bjtmp = bj + bi[i]; 590c4413a7SShri Abhyankar nzL = bi[i+1] - bi[i]; 600c4413a7SShri Abhyankar for (k=0; k < nzL; k++) { 610c4413a7SShri Abhyankar row = bjtmp[k]; 620c4413a7SShri Abhyankar pc = rtmp + bs2*row; 63c35f09e5SBarry Smith for (flg=0,j=0; j<bs2; j++) { 64c35f09e5SBarry Smith if (pc[j] != (PetscScalar)0.0) { 65c35f09e5SBarry Smith flg = 1; 66c35f09e5SBarry Smith break; 67c35f09e5SBarry Smith } 68c35f09e5SBarry Smith } 690c4413a7SShri Abhyankar if (flg) { 700c4413a7SShri Abhyankar pv = b->a + bs2*bdiag[row]; 7196b95a6bSBarry Smith /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 7296b95a6bSBarry Smith ierr = PetscKernel_A_gets_A_times_B_2(pc,pv,mwork);CHKERRQ(ierr); 730c4413a7SShri Abhyankar 740c4413a7SShri Abhyankar pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */ 750c4413a7SShri Abhyankar pv = b->a + bs2*(bdiag[row+1]+1); 760c4413a7SShri Abhyankar nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */ 770c4413a7SShri Abhyankar for (j=0; j<nz; j++) { 7896b95a6bSBarry Smith /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 790c4413a7SShri Abhyankar /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 800c4413a7SShri Abhyankar v = rtmp + 4*pj[j]; 8196b95a6bSBarry Smith ierr = PetscKernel_A_gets_A_minus_B_times_C_2(v,pc,pv);CHKERRQ(ierr); 820c4413a7SShri Abhyankar pv += 4; 830c4413a7SShri Abhyankar } 840c4413a7SShri Abhyankar ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 850c4413a7SShri Abhyankar } 860c4413a7SShri Abhyankar } 870c4413a7SShri Abhyankar 880c4413a7SShri Abhyankar /* finished row so stick it into b->a */ 890c4413a7SShri Abhyankar /* L part */ 900c4413a7SShri Abhyankar pv = b->a + bs2*bi[i]; 910c4413a7SShri Abhyankar pj = b->j + bi[i]; 920c4413a7SShri Abhyankar nz = bi[i+1] - bi[i]; 930c4413a7SShri Abhyankar for (j=0; j<nz; j++) { 940c4413a7SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 950c4413a7SShri Abhyankar } 960c4413a7SShri Abhyankar 970c4413a7SShri Abhyankar /* Mark diagonal and invert diagonal for simplier triangular solves */ 980c4413a7SShri Abhyankar pv = b->a + bs2*bdiag[i]; 990c4413a7SShri Abhyankar pj = b->j + bdiag[i]; 1000c4413a7SShri Abhyankar ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 10196b95a6bSBarry Smith /* ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 10296b95a6bSBarry Smith ierr = PetscKernel_A_gets_inverse_A_2(pv,shift);CHKERRQ(ierr); 1030c4413a7SShri Abhyankar 1040c4413a7SShri Abhyankar /* U part */ 1050c4413a7SShri Abhyankar pv = b->a + bs2*(bdiag[i+1]+1); 1060c4413a7SShri Abhyankar pj = b->j + bdiag[i+1]+1; 1070c4413a7SShri Abhyankar nz = bdiag[i] - bdiag[i+1] - 1; 1080c4413a7SShri Abhyankar for (j=0; j<nz; j++) { 1090c4413a7SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1100c4413a7SShri Abhyankar } 1110c4413a7SShri Abhyankar } 1120c4413a7SShri Abhyankar 113fca92195SBarry Smith ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); 1140c4413a7SShri Abhyankar ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 1150c4413a7SShri Abhyankar ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 11626fbe8dcSKarl Rupp 1174dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_2; 1184dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2; 1190c4413a7SShri Abhyankar C->assembled = PETSC_TRUE; 12026fbe8dcSKarl Rupp 121766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*2*2*2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 1220c4413a7SShri Abhyankar PetscFunctionReturn(0); 1230c4413a7SShri Abhyankar } 1240c4413a7SShri Abhyankar 1254c000e2eSHong Zhang #undef __FUNCT__ 1264dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering" 1274dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info) 1284c000e2eSHong Zhang { 1294c000e2eSHong Zhang Mat C =B; 1304c000e2eSHong Zhang Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data; 1314c000e2eSHong Zhang PetscErrorCode ierr; 132bbd65245SShri Abhyankar PetscInt i,j,k,nz,nzL,row,*pj; 133bbd65245SShri Abhyankar const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,bs2=a->bs2; 134bbd65245SShri Abhyankar const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag; 135bbd65245SShri Abhyankar MatScalar *rtmp,*pc,*mwork,*pv; 136bbd65245SShri Abhyankar MatScalar *aa=a->a,*v; 137bbd65245SShri Abhyankar PetscInt flg; 138182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 1394c000e2eSHong Zhang 1404c000e2eSHong Zhang PetscFunctionBegin; 1414c000e2eSHong Zhang /* generate work space needed by the factorization */ 142dcca6d9dSJed Brown ierr = PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);CHKERRQ(ierr); 1434c000e2eSHong Zhang ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 1444c000e2eSHong Zhang 1454c000e2eSHong Zhang for (i=0; i<n; i++) { 1464c000e2eSHong Zhang /* zero rtmp */ 1474c000e2eSHong Zhang /* L part */ 1484c000e2eSHong Zhang nz = bi[i+1] - bi[i]; 1494c000e2eSHong Zhang bjtmp = bj + bi[i]; 1504c000e2eSHong Zhang for (j=0; j<nz; j++) { 1514c000e2eSHong Zhang ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1524c000e2eSHong Zhang } 1534c000e2eSHong Zhang 1544c000e2eSHong Zhang /* U part */ 155b2b2dd24SShri Abhyankar nz = bdiag[i] - bdiag[i+1]; 156b2b2dd24SShri Abhyankar bjtmp = bj + bdiag[i+1]+1; 157b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 158b2b2dd24SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 159b2b2dd24SShri Abhyankar } 160b2b2dd24SShri Abhyankar 161b2b2dd24SShri Abhyankar /* load in initial (unfactored row) */ 162b2b2dd24SShri Abhyankar nz = ai[i+1] - ai[i]; 163b2b2dd24SShri Abhyankar ajtmp = aj + ai[i]; 164b2b2dd24SShri Abhyankar v = aa + bs2*ai[i]; 165b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 166b2b2dd24SShri Abhyankar ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 167b2b2dd24SShri Abhyankar } 168b2b2dd24SShri Abhyankar 169b2b2dd24SShri Abhyankar /* elimination */ 170b2b2dd24SShri Abhyankar bjtmp = bj + bi[i]; 171b2b2dd24SShri Abhyankar nzL = bi[i+1] - bi[i]; 172b2b2dd24SShri Abhyankar for (k=0; k < nzL; k++) { 173b2b2dd24SShri Abhyankar row = bjtmp[k]; 174b2b2dd24SShri Abhyankar pc = rtmp + bs2*row; 175c35f09e5SBarry Smith for (flg=0,j=0; j<bs2; j++) { 176c35f09e5SBarry Smith if (pc[j]!=(PetscScalar)0.0) { 177c35f09e5SBarry Smith flg = 1; 178c35f09e5SBarry Smith break; 179c35f09e5SBarry Smith } 180c35f09e5SBarry Smith } 181b2b2dd24SShri Abhyankar if (flg) { 182b2b2dd24SShri Abhyankar pv = b->a + bs2*bdiag[row]; 18396b95a6bSBarry Smith /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 18496b95a6bSBarry Smith ierr = PetscKernel_A_gets_A_times_B_2(pc,pv,mwork);CHKERRQ(ierr); 185b2b2dd24SShri Abhyankar 186b2b2dd24SShri Abhyankar pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */ 187b2b2dd24SShri Abhyankar pv = b->a + bs2*(bdiag[row+1]+1); 188b2b2dd24SShri Abhyankar nz = bdiag[row]-bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */ 189b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 19096b95a6bSBarry Smith /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 191b2b2dd24SShri Abhyankar /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 192b2b2dd24SShri Abhyankar v = rtmp + 4*pj[j]; 19396b95a6bSBarry Smith ierr = PetscKernel_A_gets_A_minus_B_times_C_2(v,pc,pv);CHKERRQ(ierr); 194b2b2dd24SShri Abhyankar pv += 4; 195b2b2dd24SShri Abhyankar } 196b2b2dd24SShri Abhyankar ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 197b2b2dd24SShri Abhyankar } 198b2b2dd24SShri Abhyankar } 199b2b2dd24SShri Abhyankar 200b2b2dd24SShri Abhyankar /* finished row so stick it into b->a */ 201b2b2dd24SShri Abhyankar /* L part */ 202b2b2dd24SShri Abhyankar pv = b->a + bs2*bi[i]; 203b2b2dd24SShri Abhyankar pj = b->j + bi[i]; 204b2b2dd24SShri Abhyankar nz = bi[i+1] - bi[i]; 205b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 206b2b2dd24SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 207b2b2dd24SShri Abhyankar } 208b2b2dd24SShri Abhyankar 209b2b2dd24SShri Abhyankar /* Mark diagonal and invert diagonal for simplier triangular solves */ 210b2b2dd24SShri Abhyankar pv = b->a + bs2*bdiag[i]; 211b2b2dd24SShri Abhyankar pj = b->j + bdiag[i]; 212b2b2dd24SShri Abhyankar ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 21396b95a6bSBarry Smith /* ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 21496b95a6bSBarry Smith ierr = PetscKernel_A_gets_inverse_A_2(pv,shift);CHKERRQ(ierr); 215b2b2dd24SShri Abhyankar 216b2b2dd24SShri Abhyankar /* U part */ 217b2b2dd24SShri Abhyankar /* 218b2b2dd24SShri Abhyankar pv = b->a + bs2*bi[2*n-i]; 219b2b2dd24SShri Abhyankar pj = b->j + bi[2*n-i]; 220b2b2dd24SShri Abhyankar nz = bi[2*n-i+1] - bi[2*n-i] - 1; 221b2b2dd24SShri Abhyankar */ 222b2b2dd24SShri Abhyankar pv = b->a + bs2*(bdiag[i+1]+1); 223b2b2dd24SShri Abhyankar pj = b->j + bdiag[i+1]+1; 224b2b2dd24SShri Abhyankar nz = bdiag[i] - bdiag[i+1] - 1; 225b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 226b2b2dd24SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 227b2b2dd24SShri Abhyankar } 228b2b2dd24SShri Abhyankar } 229fca92195SBarry Smith ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); 230ae3d28f0SHong Zhang 2314dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering; 232*9f5c0bcdSBarry Smith C->ops->forwardsolve = MatForwardSolve_SeqBAIJ_2_NaturalOrdering; 233*9f5c0bcdSBarry Smith C->ops->backwardsolve = MatBackwardSolve_SeqBAIJ_2_NaturalOrdering; 2344dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering; 235b2b2dd24SShri Abhyankar C->assembled = PETSC_TRUE; 23626fbe8dcSKarl Rupp 237766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*2*2*2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 238b2b2dd24SShri Abhyankar PetscFunctionReturn(0); 239b2b2dd24SShri Abhyankar } 240b2b2dd24SShri Abhyankar 241b2b2dd24SShri Abhyankar #undef __FUNCT__ 24206e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_inplace" 24306e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_inplace(Mat B,Mat A,const MatFactorInfo *info) 2444eeb42bcSBarry Smith { 245719d5645SBarry Smith Mat C = B; 2464eeb42bcSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data; 2477cdf2dbbSSatish Balay IS isrow = b->row,isicol = b->icol; 2486849ba73SBarry Smith PetscErrorCode ierr; 2495d0c19d7SBarry Smith const PetscInt *r,*ic; 2505d0c19d7SBarry Smith PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 251690b6cddSBarry Smith PetscInt *ajtmpold,*ajtmp,nz,row; 252690b6cddSBarry Smith PetscInt *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj; 253329f5518SBarry Smith MatScalar *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4; 2542515f8d2SSatish Balay MatScalar p1,p2,p3,p4; 2553f1db9ecSBarry Smith MatScalar *ba = b->a,*aa = a->a; 256182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 2574eeb42bcSBarry Smith 2583a40ed3dSBarry Smith PetscFunctionBegin; 2594eeb42bcSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 2604eeb42bcSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 261785e854fSJed Brown ierr = PetscMalloc1(4*(n+1),&rtmp);CHKERRQ(ierr); 2624eeb42bcSBarry Smith 2634eeb42bcSBarry Smith for (i=0; i<n; i++) { 2644078e994SBarry Smith nz = bi[i+1] - bi[i]; 2654078e994SBarry Smith ajtmp = bj + bi[i]; 2664eeb42bcSBarry Smith for (j=0; j<nz; j++) { 2674eeb42bcSBarry Smith x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0; 2684eeb42bcSBarry Smith } 2694eeb42bcSBarry Smith /* load in initial (unfactored row) */ 2704eeb42bcSBarry Smith idx = r[i]; 2714078e994SBarry Smith nz = ai[idx+1] - ai[idx]; 2724078e994SBarry Smith ajtmpold = aj + ai[idx]; 2734078e994SBarry Smith v = aa + 4*ai[idx]; 2744eeb42bcSBarry Smith for (j=0; j<nz; j++) { 2754eeb42bcSBarry Smith x = rtmp+4*ic[ajtmpold[j]]; 2764eeb42bcSBarry Smith x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 2774eeb42bcSBarry Smith v += 4; 2784eeb42bcSBarry Smith } 2794eeb42bcSBarry Smith row = *ajtmp++; 2804eeb42bcSBarry Smith while (row < i) { 2814eeb42bcSBarry Smith pc = rtmp + 4*row; 2824eeb42bcSBarry Smith p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 283d4a378daSJed Brown if (p1 != (PetscScalar)0.0 || p2 != (PetscScalar)0.0 || p3 != (PetscScalar)0.0 || p4 != (PetscScalar)0.0) { 2844078e994SBarry Smith pv = ba + 4*diag_offset[row]; 2854078e994SBarry Smith pj = bj + diag_offset[row] + 1; 2864eeb42bcSBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 2874eeb42bcSBarry Smith pc[0] = m1 = p1*x1 + p3*x2; 2884eeb42bcSBarry Smith pc[1] = m2 = p2*x1 + p4*x2; 2894eeb42bcSBarry Smith pc[2] = m3 = p1*x3 + p3*x4; 2904eeb42bcSBarry Smith pc[3] = m4 = p2*x3 + p4*x4; 2914078e994SBarry Smith nz = bi[row+1] - diag_offset[row] - 1; 2924eeb42bcSBarry Smith pv += 4; 2934eeb42bcSBarry Smith for (j=0; j<nz; j++) { 2944eeb42bcSBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 2954eeb42bcSBarry Smith x = rtmp + 4*pj[j]; 2964eeb42bcSBarry Smith x[0] -= m1*x1 + m3*x2; 2974eeb42bcSBarry Smith x[1] -= m2*x1 + m4*x2; 2984eeb42bcSBarry Smith x[2] -= m1*x3 + m3*x4; 2994eeb42bcSBarry Smith x[3] -= m2*x3 + m4*x4; 3004eeb42bcSBarry Smith pv += 4; 3014eeb42bcSBarry Smith } 302dc0b31edSSatish Balay ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr); 3034eeb42bcSBarry Smith } 3044eeb42bcSBarry Smith row = *ajtmp++; 3054eeb42bcSBarry Smith } 3064eeb42bcSBarry Smith /* finished row so stick it into b->a */ 3074078e994SBarry Smith pv = ba + 4*bi[i]; 3084078e994SBarry Smith pj = bj + bi[i]; 3094078e994SBarry Smith nz = bi[i+1] - bi[i]; 3104eeb42bcSBarry Smith for (j=0; j<nz; j++) { 3114eeb42bcSBarry Smith x = rtmp+4*pj[j]; 3124eeb42bcSBarry Smith pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 3134eeb42bcSBarry Smith pv += 4; 3144eeb42bcSBarry Smith } 3154eeb42bcSBarry Smith /* invert diagonal block */ 3164078e994SBarry Smith w = ba + 4*diag_offset[i]; 31796b95a6bSBarry Smith ierr = PetscKernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr); 3184eeb42bcSBarry Smith } 3194eeb42bcSBarry Smith 320606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 3214eeb42bcSBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 3224eeb42bcSBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 32326fbe8dcSKarl Rupp 32406e38f1dSHong Zhang C->ops->solve = MatSolve_SeqBAIJ_2_inplace; 32506e38f1dSHong Zhang C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_inplace; 3264eeb42bcSBarry Smith C->assembled = PETSC_TRUE; 32726fbe8dcSKarl Rupp 328766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 3293a40ed3dSBarry Smith PetscFunctionReturn(0); 3304eeb42bcSBarry Smith } 3314cd38560SBarry Smith /* 3324cd38560SBarry Smith Version for when blocks are 2 by 2 Using natural ordering 3334cd38560SBarry Smith */ 3344a2ae208SSatish Balay #undef __FUNCT__ 33506e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace" 33606e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info) 3374cd38560SBarry Smith { 3384cd38560SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data; 339dfbe8321SBarry Smith PetscErrorCode ierr; 340690b6cddSBarry Smith PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 341690b6cddSBarry Smith PetscInt *ajtmpold,*ajtmp,nz,row; 342690b6cddSBarry Smith PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; 343329f5518SBarry Smith MatScalar *pv,*v,*rtmp,*pc,*w,*x; 3442515f8d2SSatish Balay MatScalar p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4; 3454cd38560SBarry Smith MatScalar *ba = b->a,*aa = a->a; 346182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 3474cd38560SBarry Smith 3484cd38560SBarry Smith PetscFunctionBegin; 349785e854fSJed Brown ierr = PetscMalloc1(4*(n+1),&rtmp);CHKERRQ(ierr); 3504cd38560SBarry Smith for (i=0; i<n; i++) { 3514cd38560SBarry Smith nz = bi[i+1] - bi[i]; 3524cd38560SBarry Smith ajtmp = bj + bi[i]; 3534cd38560SBarry Smith for (j=0; j<nz; j++) { 3544cd38560SBarry Smith x = rtmp+4*ajtmp[j]; 3554cd38560SBarry Smith x[0] = x[1] = x[2] = x[3] = 0.0; 3564cd38560SBarry Smith } 3574cd38560SBarry Smith /* load in initial (unfactored row) */ 3584cd38560SBarry Smith nz = ai[i+1] - ai[i]; 3594cd38560SBarry Smith ajtmpold = aj + ai[i]; 3604cd38560SBarry Smith v = aa + 4*ai[i]; 3614cd38560SBarry Smith for (j=0; j<nz; j++) { 3624cd38560SBarry Smith x = rtmp+4*ajtmpold[j]; 3634cd38560SBarry Smith x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 3644cd38560SBarry Smith v += 4; 3654cd38560SBarry Smith } 3664cd38560SBarry Smith row = *ajtmp++; 3674cd38560SBarry Smith while (row < i) { 3684cd38560SBarry Smith pc = rtmp + 4*row; 3694cd38560SBarry Smith p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 370d4a378daSJed Brown if (p1 != (PetscScalar)0.0 || p2 != (PetscScalar)0.0 || p3 != (PetscScalar)0.0 || p4 != (PetscScalar)0.0) { 3714cd38560SBarry Smith pv = ba + 4*diag_offset[row]; 3724cd38560SBarry Smith pj = bj + diag_offset[row] + 1; 3734cd38560SBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 3744cd38560SBarry Smith pc[0] = m1 = p1*x1 + p3*x2; 3754cd38560SBarry Smith pc[1] = m2 = p2*x1 + p4*x2; 3764cd38560SBarry Smith pc[2] = m3 = p1*x3 + p3*x4; 3774cd38560SBarry Smith pc[3] = m4 = p2*x3 + p4*x4; 3784cd38560SBarry Smith nz = bi[row+1] - diag_offset[row] - 1; 3794cd38560SBarry Smith pv += 4; 3804cd38560SBarry Smith for (j=0; j<nz; j++) { 3814cd38560SBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 3824cd38560SBarry Smith x = rtmp + 4*pj[j]; 3834cd38560SBarry Smith x[0] -= m1*x1 + m3*x2; 3844cd38560SBarry Smith x[1] -= m2*x1 + m4*x2; 3854cd38560SBarry Smith x[2] -= m1*x3 + m3*x4; 3864cd38560SBarry Smith x[3] -= m2*x3 + m4*x4; 3874cd38560SBarry Smith pv += 4; 3884cd38560SBarry Smith } 389dc0b31edSSatish Balay ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr); 3904cd38560SBarry Smith } 3914cd38560SBarry Smith row = *ajtmp++; 3924cd38560SBarry Smith } 3934cd38560SBarry Smith /* finished row so stick it into b->a */ 3944cd38560SBarry Smith pv = ba + 4*bi[i]; 3954cd38560SBarry Smith pj = bj + bi[i]; 3964cd38560SBarry Smith nz = bi[i+1] - bi[i]; 3974cd38560SBarry Smith for (j=0; j<nz; j++) { 3984cd38560SBarry Smith x = rtmp+4*pj[j]; 3994cd38560SBarry Smith pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 4002efa7f71SHong Zhang /* 4012efa7f71SHong Zhang printf(" col %d:",pj[j]); 4022efa7f71SHong Zhang PetscInt j1; 4032efa7f71SHong Zhang for (j1=0; j1<4; j1++) printf(" %g,",*(pv+j1)); 4042efa7f71SHong Zhang printf("\n"); 4052efa7f71SHong Zhang */ 4064cd38560SBarry Smith pv += 4; 4074cd38560SBarry Smith } 4084cd38560SBarry Smith /* invert diagonal block */ 4094cd38560SBarry Smith w = ba + 4*diag_offset[i]; 4102efa7f71SHong Zhang /* 4112efa7f71SHong Zhang printf(" \n%d -th: diag: ",i); 4122efa7f71SHong Zhang for (j=0; j<4; j++) { 4132efa7f71SHong Zhang printf(" %g,",w[j]); 4142efa7f71SHong Zhang } 4152efa7f71SHong Zhang printf("\n----------------------------\n"); 4162efa7f71SHong Zhang */ 41796b95a6bSBarry Smith ierr = PetscKernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr); 4184cd38560SBarry Smith } 4194cd38560SBarry Smith 420606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 42126fbe8dcSKarl Rupp 42206e38f1dSHong Zhang C->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering_inplace; 42306e38f1dSHong Zhang C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering_inplace; 4244cd38560SBarry Smith C->assembled = PETSC_TRUE; 42526fbe8dcSKarl Rupp 426766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 4274cd38560SBarry Smith PetscFunctionReturn(0); 4284cd38560SBarry Smith } 4297fc0212eSBarry Smith 4307fc0212eSBarry Smith /* ----------------------------------------------------------- */ 4317fc0212eSBarry Smith /* 4327fc0212eSBarry Smith Version for when blocks are 1 by 1. 4337fc0212eSBarry Smith */ 4344a2ae208SSatish Balay #undef __FUNCT__ 435048b5e81SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1" 436048b5e81SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat B,Mat A,const MatFactorInfo *info) 437048b5e81SShri Abhyankar { 438048b5e81SShri Abhyankar Mat C =B; 439048b5e81SShri Abhyankar Mat_SeqBAIJ *a =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data; 440048b5e81SShri Abhyankar IS isrow = b->row,isicol = b->icol; 441048b5e81SShri Abhyankar PetscErrorCode ierr; 442048b5e81SShri Abhyankar const PetscInt *r,*ic,*ics; 443048b5e81SShri Abhyankar const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bdiag=b->diag; 444048b5e81SShri Abhyankar PetscInt i,j,k,nz,nzL,row,*pj; 445048b5e81SShri Abhyankar const PetscInt *ajtmp,*bjtmp; 446048b5e81SShri Abhyankar MatScalar *rtmp,*pc,multiplier,*pv; 447048b5e81SShri Abhyankar const MatScalar *aa=a->a,*v; 448ace3abfcSBarry Smith PetscBool row_identity,col_identity; 449048b5e81SShri Abhyankar FactorShiftCtx sctx; 450048b5e81SShri Abhyankar const PetscInt *ddiag; 451048b5e81SShri Abhyankar PetscReal rs; 452048b5e81SShri Abhyankar MatScalar d; 453048b5e81SShri Abhyankar 454048b5e81SShri Abhyankar PetscFunctionBegin; 455048b5e81SShri Abhyankar /* MatPivotSetUp(): initialize shift context sctx */ 456048b5e81SShri Abhyankar ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr); 457048b5e81SShri Abhyankar 458048b5e81SShri Abhyankar if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */ 459048b5e81SShri Abhyankar ddiag = a->diag; 460048b5e81SShri Abhyankar sctx.shift_top = info->zeropivot; 461048b5e81SShri Abhyankar for (i=0; i<n; i++) { 462048b5e81SShri Abhyankar /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */ 463048b5e81SShri Abhyankar d = (aa)[ddiag[i]]; 464048b5e81SShri Abhyankar rs = -PetscAbsScalar(d) - PetscRealPart(d); 465048b5e81SShri Abhyankar v = aa+ai[i]; 466048b5e81SShri Abhyankar nz = ai[i+1] - ai[i]; 46726fbe8dcSKarl Rupp for (j=0; j<nz; j++) rs += PetscAbsScalar(v[j]); 468048b5e81SShri Abhyankar if (rs>sctx.shift_top) sctx.shift_top = rs; 469048b5e81SShri Abhyankar } 470048b5e81SShri Abhyankar sctx.shift_top *= 1.1; 471048b5e81SShri Abhyankar sctx.nshift_max = 5; 472048b5e81SShri Abhyankar sctx.shift_lo = 0.; 473048b5e81SShri Abhyankar sctx.shift_hi = 1.; 474048b5e81SShri Abhyankar } 475048b5e81SShri Abhyankar 476048b5e81SShri Abhyankar ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 477048b5e81SShri Abhyankar ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 478785e854fSJed Brown ierr = PetscMalloc1((n+1),&rtmp);CHKERRQ(ierr); 479048b5e81SShri Abhyankar ics = ic; 480048b5e81SShri Abhyankar 481048b5e81SShri Abhyankar do { 482048b5e81SShri Abhyankar sctx.newshift = PETSC_FALSE; 483048b5e81SShri Abhyankar for (i=0; i<n; i++) { 484048b5e81SShri Abhyankar /* zero rtmp */ 485048b5e81SShri Abhyankar /* L part */ 486048b5e81SShri Abhyankar nz = bi[i+1] - bi[i]; 487048b5e81SShri Abhyankar bjtmp = bj + bi[i]; 488048b5e81SShri Abhyankar for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0; 489048b5e81SShri Abhyankar 490048b5e81SShri Abhyankar /* U part */ 491048b5e81SShri Abhyankar nz = bdiag[i]-bdiag[i+1]; 492048b5e81SShri Abhyankar bjtmp = bj + bdiag[i+1]+1; 493048b5e81SShri Abhyankar for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0; 494048b5e81SShri Abhyankar 495048b5e81SShri Abhyankar /* load in initial (unfactored row) */ 496048b5e81SShri Abhyankar nz = ai[r[i]+1] - ai[r[i]]; 497048b5e81SShri Abhyankar ajtmp = aj + ai[r[i]]; 498048b5e81SShri Abhyankar v = aa + ai[r[i]]; 49926fbe8dcSKarl Rupp for (j=0; j<nz; j++) rtmp[ics[ajtmp[j]]] = v[j]; 50026fbe8dcSKarl Rupp 501048b5e81SShri Abhyankar /* ZeropivotApply() */ 502048b5e81SShri Abhyankar rtmp[i] += sctx.shift_amount; /* shift the diagonal of the matrix */ 503048b5e81SShri Abhyankar 504048b5e81SShri Abhyankar /* elimination */ 505048b5e81SShri Abhyankar bjtmp = bj + bi[i]; 506048b5e81SShri Abhyankar row = *bjtmp++; 507048b5e81SShri Abhyankar nzL = bi[i+1] - bi[i]; 508048b5e81SShri Abhyankar for (k=0; k < nzL; k++) { 509048b5e81SShri Abhyankar pc = rtmp + row; 510d4a378daSJed Brown if (*pc != (PetscScalar)0.0) { 511048b5e81SShri Abhyankar pv = b->a + bdiag[row]; 512048b5e81SShri Abhyankar multiplier = *pc * (*pv); 513048b5e81SShri Abhyankar *pc = multiplier; 51426fbe8dcSKarl Rupp 515048b5e81SShri Abhyankar pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */ 516048b5e81SShri Abhyankar pv = b->a + bdiag[row+1]+1; 517048b5e81SShri Abhyankar nz = bdiag[row]-bdiag[row+1]-1; /* num of entries in U(row,:) excluding diag */ 518048b5e81SShri Abhyankar for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j]; 519048b5e81SShri Abhyankar ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); 520048b5e81SShri Abhyankar } 521048b5e81SShri Abhyankar row = *bjtmp++; 522048b5e81SShri Abhyankar } 523048b5e81SShri Abhyankar 524048b5e81SShri Abhyankar /* finished row so stick it into b->a */ 525048b5e81SShri Abhyankar rs = 0.0; 526048b5e81SShri Abhyankar /* L part */ 527048b5e81SShri Abhyankar pv = b->a + bi[i]; 528048b5e81SShri Abhyankar pj = b->j + bi[i]; 529048b5e81SShri Abhyankar nz = bi[i+1] - bi[i]; 530048b5e81SShri Abhyankar for (j=0; j<nz; j++) { 531048b5e81SShri Abhyankar pv[j] = rtmp[pj[j]]; rs += PetscAbsScalar(pv[j]); 532048b5e81SShri Abhyankar } 533048b5e81SShri Abhyankar 534048b5e81SShri Abhyankar /* U part */ 535048b5e81SShri Abhyankar pv = b->a + bdiag[i+1]+1; 536048b5e81SShri Abhyankar pj = b->j + bdiag[i+1]+1; 537048b5e81SShri Abhyankar nz = bdiag[i] - bdiag[i+1]-1; 538048b5e81SShri Abhyankar for (j=0; j<nz; j++) { 539048b5e81SShri Abhyankar pv[j] = rtmp[pj[j]]; rs += PetscAbsScalar(pv[j]); 540048b5e81SShri Abhyankar } 541048b5e81SShri Abhyankar 542048b5e81SShri Abhyankar sctx.rs = rs; 543048b5e81SShri Abhyankar sctx.pv = rtmp[i]; 544d0660788SBarry Smith ierr = MatPivotCheck(A,info,&sctx,i);CHKERRQ(ierr); 545048b5e81SShri Abhyankar if (sctx.newshift) break; /* break for-loop */ 546048b5e81SShri Abhyankar rtmp[i] = sctx.pv; /* sctx.pv might be updated in the case of MAT_SHIFT_INBLOCKS */ 547048b5e81SShri Abhyankar 548048b5e81SShri Abhyankar /* Mark diagonal and invert diagonal for simplier triangular solves */ 549048b5e81SShri Abhyankar pv = b->a + bdiag[i]; 550d4a378daSJed Brown *pv = (PetscScalar)1.0/rtmp[i]; 551048b5e81SShri Abhyankar 552048b5e81SShri Abhyankar } /* endof for (i=0; i<n; i++) { */ 553048b5e81SShri Abhyankar 554048b5e81SShri Abhyankar /* MatPivotRefine() */ 555048b5e81SShri Abhyankar if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE && !sctx.newshift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) { 556048b5e81SShri Abhyankar /* 557048b5e81SShri Abhyankar * if no shift in this attempt & shifting & started shifting & can refine, 558048b5e81SShri Abhyankar * then try lower shift 559048b5e81SShri Abhyankar */ 560048b5e81SShri Abhyankar sctx.shift_hi = sctx.shift_fraction; 561048b5e81SShri Abhyankar sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.; 562048b5e81SShri Abhyankar sctx.shift_amount = sctx.shift_fraction * sctx.shift_top; 563048b5e81SShri Abhyankar sctx.newshift = PETSC_TRUE; 564048b5e81SShri Abhyankar sctx.nshift++; 565048b5e81SShri Abhyankar } 566048b5e81SShri Abhyankar } while (sctx.newshift); 567048b5e81SShri Abhyankar 568048b5e81SShri Abhyankar ierr = PetscFree(rtmp);CHKERRQ(ierr); 569048b5e81SShri Abhyankar ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 570048b5e81SShri Abhyankar ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 571048b5e81SShri Abhyankar 572048b5e81SShri Abhyankar ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 573048b5e81SShri Abhyankar ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr); 574048b5e81SShri Abhyankar if (row_identity && col_identity) { 575048b5e81SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_1_NaturalOrdering; 576*9f5c0bcdSBarry Smith C->ops->forwardsolve = MatForwardSolve_SeqBAIJ_1_NaturalOrdering; 577*9f5c0bcdSBarry Smith C->ops->backwardsolve = MatBackwardSolve_SeqBAIJ_1_NaturalOrdering; 57893fd935bSShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering; 579048b5e81SShri Abhyankar } else { 580048b5e81SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_1; 58193fd935bSShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1; 582048b5e81SShri Abhyankar } 583048b5e81SShri Abhyankar C->assembled = PETSC_TRUE; 584048b5e81SShri Abhyankar ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr); 585048b5e81SShri Abhyankar 586048b5e81SShri Abhyankar /* MatShiftView(A,info,&sctx) */ 587048b5e81SShri Abhyankar if (sctx.nshift) { 588048b5e81SShri Abhyankar if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { 58957622a8eSBarry Smith ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %g, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,(double)sctx.shift_amount,(double)sctx.shift_fraction,(double)sctx.shift_top);CHKERRQ(ierr); 590048b5e81SShri Abhyankar } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 59157622a8eSBarry Smith ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr); 592048b5e81SShri Abhyankar } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) { 59357622a8eSBarry Smith ierr = PetscInfo2(A,"number of shift_inblocks applied %D, each shift_amount %g\n",sctx.nshift,(double)info->shiftamount);CHKERRQ(ierr); 594048b5e81SShri Abhyankar } 595048b5e81SShri Abhyankar } 596048b5e81SShri Abhyankar PetscFunctionReturn(0); 597048b5e81SShri Abhyankar } 598048b5e81SShri Abhyankar 599048b5e81SShri Abhyankar #undef __FUNCT__ 60006e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1_inplace" 60106e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info) 6027fc0212eSBarry Smith { 6037fc0212eSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data; 6047cdf2dbbSSatish Balay IS isrow = b->row,isicol = b->icol; 6056849ba73SBarry Smith PetscErrorCode ierr; 6065d0c19d7SBarry Smith const PetscInt *r,*ic; 6075d0c19d7SBarry Smith PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 608690b6cddSBarry Smith PetscInt *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j; 609690b6cddSBarry Smith PetscInt *diag_offset = b->diag,diag,*pj; 610329f5518SBarry Smith MatScalar *pv,*v,*rtmp,multiplier,*pc; 6113f1db9ecSBarry Smith MatScalar *ba = b->a,*aa = a->a; 612ace3abfcSBarry Smith PetscBool row_identity, col_identity; 6137fc0212eSBarry Smith 6143a40ed3dSBarry Smith PetscFunctionBegin; 6157fc0212eSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 6167fc0212eSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 617785e854fSJed Brown ierr = PetscMalloc1((n+1),&rtmp);CHKERRQ(ierr); 6187fc0212eSBarry Smith 6197fc0212eSBarry Smith for (i=0; i<n; i++) { 6204078e994SBarry Smith nz = bi[i+1] - bi[i]; 6214078e994SBarry Smith ajtmp = bj + bi[i]; 6227fc0212eSBarry Smith for (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0; 6237fc0212eSBarry Smith 6247fc0212eSBarry Smith /* load in initial (unfactored row) */ 6254078e994SBarry Smith nz = ai[r[i]+1] - ai[r[i]]; 6264078e994SBarry Smith ajtmpold = aj + ai[r[i]]; 6274078e994SBarry Smith v = aa + ai[r[i]]; 6287fc0212eSBarry Smith for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] = v[j]; 6297fc0212eSBarry Smith 6307fc0212eSBarry Smith row = *ajtmp++; 6317fc0212eSBarry Smith while (row < i) { 6327fc0212eSBarry Smith pc = rtmp + row; 6337fc0212eSBarry Smith if (*pc != 0.0) { 6344078e994SBarry Smith pv = ba + diag_offset[row]; 6354078e994SBarry Smith pj = bj + diag_offset[row] + 1; 6367fc0212eSBarry Smith multiplier = *pc * *pv++; 6377fc0212eSBarry Smith *pc = multiplier; 6384078e994SBarry Smith nz = bi[row+1] - diag_offset[row] - 1; 6397fc0212eSBarry Smith for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j]; 640dc0b31edSSatish Balay ierr = PetscLogFlops(1.0+2.0*nz);CHKERRQ(ierr); 6417fc0212eSBarry Smith } 6427fc0212eSBarry Smith row = *ajtmp++; 6437fc0212eSBarry Smith } 6447fc0212eSBarry Smith /* finished row so stick it into b->a */ 6454078e994SBarry Smith pv = ba + bi[i]; 6464078e994SBarry Smith pj = bj + bi[i]; 6474078e994SBarry Smith nz = bi[i+1] - bi[i]; 64826fbe8dcSKarl Rupp for (j=0; j<nz; j++) pv[j] = rtmp[pj[j]]; 6494078e994SBarry Smith diag = diag_offset[i] - bi[i]; 6507fc0212eSBarry Smith /* check pivot entry for current row */ 65165e19b50SBarry 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); 6527fc0212eSBarry Smith pv[diag] = 1.0/pv[diag]; 6537fc0212eSBarry Smith } 6547fc0212eSBarry Smith 655606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 6567fc0212eSBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 6577fc0212eSBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 658f58c8c74SBarry Smith ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 659f58c8c74SBarry Smith ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr); 660f58c8c74SBarry Smith if (row_identity && col_identity) { 66106e38f1dSHong Zhang C->ops->solve = MatSolve_SeqBAIJ_1_NaturalOrdering_inplace; 66206e38f1dSHong Zhang C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering_inplace; 663f58c8c74SBarry Smith } else { 66406e38f1dSHong Zhang C->ops->solve = MatSolve_SeqBAIJ_1_inplace; 66506e38f1dSHong Zhang C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_inplace; 666f58c8c74SBarry Smith } 6677fc0212eSBarry Smith C->assembled = PETSC_TRUE; 668d0f46423SBarry Smith ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr); 6693a40ed3dSBarry Smith PetscFunctionReturn(0); 6707fc0212eSBarry Smith } 6717fc0212eSBarry Smith 672b24902e0SBarry Smith #undef __FUNCT__ 673b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqbaij_petsc" 6748cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A,MatFactorType ftype,Mat *B) 675b24902e0SBarry Smith { 676d0f46423SBarry Smith PetscInt n = A->rmap->n; 677b24902e0SBarry Smith PetscErrorCode ierr; 678b24902e0SBarry Smith 679b24902e0SBarry Smith PetscFunctionBegin; 680b499a2aeSHong Zhang #if defined(PETSC_USE_COMPLEX) 681b499a2aeSHong Zhang if (A->hermitian && (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Hermitian Factor is not supported"); 682b499a2aeSHong Zhang #endif 683ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),B);CHKERRQ(ierr); 684b24902e0SBarry Smith ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 685c8342467SHong Zhang if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) { 686b24902e0SBarry Smith ierr = MatSetType(*B,MATSEQBAIJ);CHKERRQ(ierr); 68726fbe8dcSKarl Rupp 6884dd39f65SShri Abhyankar (*B)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqBAIJ; 6894dd39f65SShri Abhyankar (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqBAIJ; 690b24902e0SBarry Smith } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 691b24902e0SBarry Smith ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr); 6920298fd71SBarry Smith ierr = MatSeqSBAIJSetPreallocation(*B,A->rmap->bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); 69326fbe8dcSKarl Rupp 6945c9eb25fSBarry Smith (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqBAIJ; 6955c9eb25fSBarry Smith (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ; 696e32f2f54SBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 697d5f3da31SBarry Smith (*B)->factortype = ftype; 698b24902e0SBarry Smith PetscFunctionReturn(0); 699b24902e0SBarry Smith } 700273d9f13SBarry Smith 701db4efbfdSBarry Smith #undef __FUNCT__ 70253c77d0aSJed Brown #define __FUNCT__ "MatGetFactorAvailable_seqbaij_petsc" 703ace3abfcSBarry Smith PetscErrorCode MatGetFactorAvailable_seqbaij_petsc(Mat A,MatFactorType ftype,PetscBool *flg) 704db4efbfdSBarry Smith { 705db4efbfdSBarry Smith PetscFunctionBegin; 706db4efbfdSBarry Smith *flg = PETSC_TRUE; 707db4efbfdSBarry Smith PetscFunctionReturn(0); 708db4efbfdSBarry Smith } 709db4efbfdSBarry Smith 7105d517e7eSBarry Smith /* ----------------------------------------------------------- */ 7114a2ae208SSatish Balay #undef __FUNCT__ 7124a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqBAIJ" 7130481f469SBarry Smith PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,const MatFactorInfo *info) 7145d517e7eSBarry Smith { 715dfbe8321SBarry Smith PetscErrorCode ierr; 7165d517e7eSBarry Smith Mat C; 7175d517e7eSBarry Smith 7183a40ed3dSBarry Smith PetscFunctionBegin; 7192692d6eeSBarry Smith ierr = MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr); 720719d5645SBarry Smith ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr); 721719d5645SBarry Smith ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr); 72226fbe8dcSKarl Rupp 723db4efbfdSBarry Smith A->ops->solve = C->ops->solve; 724db4efbfdSBarry Smith A->ops->solvetranspose = C->ops->solvetranspose; 72526fbe8dcSKarl Rupp 726eb6b5d47SBarry Smith ierr = MatHeaderMerge(A,C);CHKERRQ(ierr); 7273bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr); 7283a40ed3dSBarry Smith PetscFunctionReturn(0); 7295d517e7eSBarry Smith } 7304cd38560SBarry Smith 731c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h> 732c05c3958SHong Zhang #undef __FUNCT__ 733c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N" 7340481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat C,Mat A,const MatFactorInfo *info) 735c05c3958SHong Zhang { 736f3086b4bSHong Zhang PetscErrorCode ierr; 73778910aadSHong Zhang Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 73878910aadSHong Zhang Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 73978910aadSHong Zhang IS ip=b->row; 7405d0c19d7SBarry Smith const PetscInt *rip; 7415d0c19d7SBarry Smith PetscInt i,j,mbs=a->mbs,bs=A->rmap->bs,*bi=b->i,*bj=b->j,*bcol; 74278910aadSHong Zhang PetscInt *ai=a->i,*aj=a->j; 74378910aadSHong Zhang PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 74478910aadSHong Zhang MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 74507b50cabSHong Zhang PetscReal rs; 7460e95ead3SHong Zhang FactorShiftCtx sctx; 74778910aadSHong Zhang 748c05c3958SHong Zhang PetscFunctionBegin; 74907b50cabSHong Zhang if (bs > 1) { /* convert A to a SBAIJ matrix and apply Cholesky factorization from it */ 7506ad2eaddSHong Zhang if (!a->sbaijMat) { 751ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 7526ad2eaddSHong Zhang } 753719d5645SBarry Smith ierr = (a->sbaijMat)->ops->choleskyfactornumeric(C,a->sbaijMat,info);CHKERRQ(ierr); 7546bf464f9SBarry Smith ierr = MatDestroy(&a->sbaijMat);CHKERRQ(ierr); 7556ad2eaddSHong Zhang PetscFunctionReturn(0); 7566ad2eaddSHong Zhang } 75778910aadSHong Zhang 75807b50cabSHong Zhang /* MatPivotSetUp(): initialize shift context sctx */ 75907b50cabSHong Zhang ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr); 7603cea4cbeSHong Zhang 7616ad2eaddSHong Zhang ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 762dcca6d9dSJed Brown ierr = PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&jl);CHKERRQ(ierr); 76378910aadSHong Zhang 76475567043SBarry Smith sctx.shift_amount = 0.; 7653cea4cbeSHong Zhang sctx.nshift = 0; 76678910aadSHong Zhang do { 76707b50cabSHong Zhang sctx.newshift = PETSC_FALSE; 76878910aadSHong Zhang for (i=0; i<mbs; i++) { 769e60cf9a0SBarry Smith rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 77078910aadSHong Zhang } 77178910aadSHong Zhang 77278910aadSHong Zhang for (k = 0; k<mbs; k++) { 77378910aadSHong Zhang bval = ba + bi[k]; 77478910aadSHong Zhang /* initialize k-th row by the perm[k]-th row of A */ 77578910aadSHong Zhang jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 77678910aadSHong Zhang for (j = jmin; j < jmax; j++) { 77778910aadSHong Zhang col = rip[aj[j]]; 77878910aadSHong Zhang if (col >= k) { /* only take upper triangular entry */ 77978910aadSHong Zhang rtmp[col] = aa[j]; 78078910aadSHong Zhang *bval++ = 0.0; /* for in-place factorization */ 78178910aadSHong Zhang } 78278910aadSHong Zhang } 7833cea4cbeSHong Zhang 7843cea4cbeSHong Zhang /* shift the diagonal of the matrix */ 7853cea4cbeSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 78678910aadSHong Zhang 78778910aadSHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 78878910aadSHong Zhang dk = rtmp[k]; 78978910aadSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 79078910aadSHong Zhang 79178910aadSHong Zhang while (i < k) { 79278910aadSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 79378910aadSHong Zhang 79478910aadSHong Zhang /* compute multiplier, update diag(k) and U(i,k) */ 79578910aadSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 79678910aadSHong Zhang uikdi = -ba[ili]*ba[bi[i]]; /* diagonal(k) */ 79778910aadSHong Zhang dk += uikdi*ba[ili]; 79878910aadSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 79978910aadSHong Zhang 80078910aadSHong Zhang /* add multiple of row i to k-th row */ 80178910aadSHong Zhang jmin = ili + 1; jmax = bi[i+1]; 80278910aadSHong Zhang if (jmin < jmax) { 80378910aadSHong Zhang for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 80478910aadSHong Zhang /* update il and jl for row i */ 80578910aadSHong Zhang il[i] = jmin; 80678910aadSHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 80778910aadSHong Zhang } 80878910aadSHong Zhang i = nexti; 80978910aadSHong Zhang } 81078910aadSHong Zhang 8113cea4cbeSHong Zhang /* shift the diagonals when zero pivot is detected */ 8123cea4cbeSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 8133cea4cbeSHong Zhang rs = 0.0; 81478910aadSHong Zhang jmin = bi[k]+1; 81578910aadSHong Zhang nz = bi[k+1] - jmin; 81678910aadSHong Zhang if (nz) { 81778910aadSHong Zhang bcol = bj + jmin; 81878910aadSHong Zhang while (nz--) { 81989f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 82089f845c8SHong Zhang bcol++; 82178910aadSHong Zhang } 82278910aadSHong Zhang } 8233cea4cbeSHong Zhang 8243cea4cbeSHong Zhang sctx.rs = rs; 8253cea4cbeSHong Zhang sctx.pv = dk; 826d0660788SBarry Smith ierr = MatPivotCheck(A,info,&sctx,k);CHKERRQ(ierr); 82707b50cabSHong Zhang if (sctx.newshift) break; 8280e95ead3SHong Zhang dk = sctx.pv; 82978910aadSHong Zhang 83078910aadSHong Zhang /* copy data into U(k,:) */ 83178910aadSHong Zhang ba[bi[k]] = 1.0/dk; /* U(k,k) */ 83278910aadSHong Zhang jmin = bi[k]+1; jmax = bi[k+1]; 83378910aadSHong Zhang if (jmin < jmax) { 83478910aadSHong Zhang for (j=jmin; j<jmax; j++) { 83578910aadSHong Zhang col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 83678910aadSHong Zhang } 83778910aadSHong Zhang /* add the k-th row into il and jl */ 83878910aadSHong Zhang il[k] = jmin; 83978910aadSHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 84078910aadSHong Zhang } 84178910aadSHong Zhang } 84207b50cabSHong Zhang } while (sctx.newshift); 843fca92195SBarry Smith ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr); 84478910aadSHong Zhang 84578910aadSHong Zhang ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 84626fbe8dcSKarl Rupp 84778910aadSHong Zhang C->assembled = PETSC_TRUE; 84878910aadSHong Zhang C->preallocated = PETSC_TRUE; 84926fbe8dcSKarl Rupp 850d0f46423SBarry Smith ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); 8513cea4cbeSHong Zhang if (sctx.nshift) { 852f4db908eSBarry Smith if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { 85357622a8eSBarry Smith ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr); 854f4db908eSBarry Smith } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 85557622a8eSBarry Smith ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr); 85678910aadSHong Zhang } 85778910aadSHong Zhang } 858c05c3958SHong Zhang PetscFunctionReturn(0); 859c05c3958SHong Zhang } 8604cd38560SBarry Smith 861c05c3958SHong Zhang #undef __FUNCT__ 862c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering" 8630481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info) 864c05c3958SHong Zhang { 86578910aadSHong Zhang Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 86678910aadSHong Zhang Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 86778910aadSHong Zhang PetscErrorCode ierr; 8683cea4cbeSHong Zhang PetscInt i,j,am=a->mbs; 86978910aadSHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 8703cea4cbeSHong Zhang PetscInt k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz; 87178910aadSHong Zhang MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval; 8720e95ead3SHong Zhang PetscReal rs; 8730e95ead3SHong Zhang FactorShiftCtx sctx; 87478910aadSHong Zhang 875c05c3958SHong Zhang PetscFunctionBegin; 8760e95ead3SHong Zhang /* MatPivotSetUp(): initialize shift context sctx */ 8770e95ead3SHong Zhang ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr); 8783cea4cbeSHong Zhang 879dcca6d9dSJed Brown ierr = PetscMalloc3(am,&rtmp,am,&il,am,&jl);CHKERRQ(ierr); 88078910aadSHong Zhang 88178910aadSHong Zhang do { 88207b50cabSHong Zhang sctx.newshift = PETSC_FALSE; 88378910aadSHong Zhang for (i=0; i<am; i++) { 884e60cf9a0SBarry Smith rtmp[i] = 0.0; jl[i] = am; il[0] = 0; 88578910aadSHong Zhang } 88678910aadSHong Zhang 88778910aadSHong Zhang for (k = 0; k<am; k++) { 88878910aadSHong Zhang /* initialize k-th row with elements nonzero in row perm(k) of A */ 88978910aadSHong Zhang nz = ai[k+1] - ai[k]; 89078910aadSHong Zhang acol = aj + ai[k]; 89178910aadSHong Zhang aval = aa + ai[k]; 89278910aadSHong Zhang bval = ba + bi[k]; 89378910aadSHong Zhang while (nz--) { 89478910aadSHong Zhang if (*acol < k) { /* skip lower triangular entries */ 89578910aadSHong Zhang acol++; aval++; 89678910aadSHong Zhang } else { 89778910aadSHong Zhang rtmp[*acol++] = *aval++; 89878910aadSHong Zhang *bval++ = 0.0; /* for in-place factorization */ 89978910aadSHong Zhang } 90078910aadSHong Zhang } 9013cea4cbeSHong Zhang 9023cea4cbeSHong Zhang /* shift the diagonal of the matrix */ 9033cea4cbeSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 90478910aadSHong Zhang 90578910aadSHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 90678910aadSHong Zhang dk = rtmp[k]; 90778910aadSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 90878910aadSHong Zhang 90978910aadSHong Zhang while (i < k) { 91078910aadSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 91178910aadSHong Zhang /* compute multiplier, update D(k) and U(i,k) */ 91278910aadSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 91378910aadSHong Zhang uikdi = -ba[ili]*ba[bi[i]]; 91478910aadSHong Zhang dk += uikdi*ba[ili]; 91578910aadSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 91678910aadSHong Zhang 91778910aadSHong Zhang /* add multiple of row i to k-th row ... */ 91878910aadSHong Zhang jmin = ili + 1; 91978910aadSHong Zhang nz = bi[i+1] - jmin; 92078910aadSHong Zhang if (nz > 0) { 92178910aadSHong Zhang bcol = bj + jmin; 92278910aadSHong Zhang bval = ba + jmin; 92378910aadSHong Zhang while (nz--) rtmp[*bcol++] += uikdi*(*bval++); 92478910aadSHong Zhang /* update il and jl for i-th row */ 92578910aadSHong Zhang il[i] = jmin; 92678910aadSHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 92778910aadSHong Zhang } 92878910aadSHong Zhang i = nexti; 92978910aadSHong Zhang } 93078910aadSHong Zhang 9313cea4cbeSHong Zhang /* shift the diagonals when zero pivot is detected */ 9323cea4cbeSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 9333cea4cbeSHong Zhang rs = 0.0; 93478910aadSHong Zhang jmin = bi[k]+1; 93578910aadSHong Zhang nz = bi[k+1] - jmin; 93678910aadSHong Zhang if (nz) { 93778910aadSHong Zhang bcol = bj + jmin; 93878910aadSHong Zhang while (nz--) { 93989f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 94089f845c8SHong Zhang bcol++; 94178910aadSHong Zhang } 94278910aadSHong Zhang } 9433cea4cbeSHong Zhang 9443cea4cbeSHong Zhang sctx.rs = rs; 9453cea4cbeSHong Zhang sctx.pv = dk; 946d0660788SBarry Smith ierr = MatPivotCheck(A,info,&sctx,k);CHKERRQ(ierr); 94707b50cabSHong Zhang if (sctx.newshift) break; /* sctx.shift_amount is updated */ 9480e95ead3SHong Zhang dk = sctx.pv; 94978910aadSHong Zhang 95078910aadSHong Zhang /* copy data into U(k,:) */ 95178910aadSHong Zhang ba[bi[k]] = 1.0/dk; 95278910aadSHong Zhang jmin = bi[k]+1; 95378910aadSHong Zhang nz = bi[k+1] - jmin; 95478910aadSHong Zhang if (nz) { 95578910aadSHong Zhang bcol = bj + jmin; 95678910aadSHong Zhang bval = ba + jmin; 95778910aadSHong Zhang while (nz--) { 95878910aadSHong Zhang *bval++ = rtmp[*bcol]; 95978910aadSHong Zhang rtmp[*bcol++] = 0.0; 96078910aadSHong Zhang } 96178910aadSHong Zhang /* add k-th row into il and jl */ 96278910aadSHong Zhang il[k] = jmin; 96378910aadSHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 96478910aadSHong Zhang } 96578910aadSHong Zhang } 96607b50cabSHong Zhang } while (sctx.newshift); 967fca92195SBarry Smith ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr); 96878910aadSHong Zhang 9690a3c351aSHong Zhang C->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 9700a3c351aSHong Zhang C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 97178910aadSHong Zhang C->assembled = PETSC_TRUE; 97278910aadSHong Zhang C->preallocated = PETSC_TRUE; 97326fbe8dcSKarl Rupp 974d0f46423SBarry Smith ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); 9753cea4cbeSHong Zhang if (sctx.nshift) { 976f4db908eSBarry Smith if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 97757622a8eSBarry Smith ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr); 978f4db908eSBarry Smith } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { 97957622a8eSBarry Smith ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr); 98078910aadSHong Zhang } 98178910aadSHong Zhang } 982c05c3958SHong Zhang PetscFunctionReturn(0); 983c05c3958SHong Zhang } 984c05c3958SHong Zhang 985c6db04a5SJed Brown #include <petscbt.h> 986c6db04a5SJed Brown #include <../src/mat/utils/freespace.h> 987c05c3958SHong Zhang #undef __FUNCT__ 988c05c3958SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ" 9890481f469SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 990c05c3958SHong Zhang { 99178910aadSHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 99278910aadSHong Zhang Mat_SeqSBAIJ *b; 99378910aadSHong Zhang Mat B; 99478910aadSHong Zhang PetscErrorCode ierr; 99548dd3d27SHong Zhang PetscBool perm_identity,missing; 9965d0c19d7SBarry Smith PetscInt reallocs=0,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap->bs,*ui; 9975d0c19d7SBarry Smith const PetscInt *rip; 99878910aadSHong Zhang PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 9990298fd71SBarry Smith PetscInt nlnk,*lnk,*lnk_lvl=NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr; 100078910aadSHong Zhang PetscReal fill =info->fill,levels=info->levels; 10010298fd71SBarry Smith PetscFreeSpaceList free_space =NULL,current_space=NULL; 10020298fd71SBarry Smith PetscFreeSpaceList free_space_lvl=NULL,current_space_lvl=NULL; 100378910aadSHong Zhang PetscBT lnkbt; 100478910aadSHong Zhang 1005c05c3958SHong Zhang PetscFunctionBegin; 100648dd3d27SHong Zhang ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr); 100748dd3d27SHong Zhang if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i); 100848dd3d27SHong Zhang 10096ad2eaddSHong Zhang if (bs > 1) { 10106ad2eaddSHong Zhang if (!a->sbaijMat) { 1011ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 10126ad2eaddSHong Zhang } 1013719d5645SBarry Smith (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */ 101426fbe8dcSKarl Rupp 1015719d5645SBarry Smith ierr = MatICCFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr); 10166ad2eaddSHong Zhang PetscFunctionReturn(0); 10176ad2eaddSHong Zhang } 10186ad2eaddSHong Zhang 101978910aadSHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 102078910aadSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 102178910aadSHong Zhang 102278910aadSHong Zhang /* special case that simply copies fill pattern */ 102378910aadSHong Zhang if (!levels && perm_identity) { 1024785e854fSJed Brown ierr = PetscMalloc1((am+1),&ui);CHKERRQ(ierr); 102526fbe8dcSKarl Rupp for (i=0; i<am; i++) ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */ 1026719d5645SBarry Smith B = fact; 102778910aadSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr); 102878910aadSHong Zhang 1029b24902e0SBarry Smith 103078910aadSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 103178910aadSHong Zhang uj = b->j; 103278910aadSHong Zhang for (i=0; i<am; i++) { 103378910aadSHong Zhang aj = a->j + a->diag[i]; 103426fbe8dcSKarl Rupp for (j=0; j<ui[i]; j++) *uj++ = *aj++; 103578910aadSHong Zhang b->ilen[i] = ui[i]; 103678910aadSHong Zhang } 103778910aadSHong Zhang ierr = PetscFree(ui);CHKERRQ(ierr); 103826fbe8dcSKarl Rupp 1039d5f3da31SBarry Smith B->factortype = MAT_FACTOR_NONE; 104026fbe8dcSKarl Rupp 104178910aadSHong Zhang ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 104278910aadSHong Zhang ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1043d5f3da31SBarry Smith B->factortype = MAT_FACTOR_ICC; 104478910aadSHong Zhang 104578910aadSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 104678910aadSHong Zhang PetscFunctionReturn(0); 104778910aadSHong Zhang } 104878910aadSHong Zhang 104978910aadSHong Zhang /* initialization */ 1050785e854fSJed Brown ierr = PetscMalloc1((am+1),&ui);CHKERRQ(ierr); 1051e60cf9a0SBarry Smith ui[0] = 0; 1052785e854fSJed Brown ierr = PetscMalloc1((2*am+1),&cols_lvl);CHKERRQ(ierr); 105378910aadSHong Zhang 105478910aadSHong Zhang /* jl: linked list for storing indices of the pivot rows 105578910aadSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 1056dcca6d9dSJed Brown ierr = PetscMalloc4(am,&uj_ptr,am,&uj_lvl_ptr,am,&il,am,&jl);CHKERRQ(ierr); 105778910aadSHong Zhang for (i=0; i<am; i++) { 1058e60cf9a0SBarry Smith jl[i] = am; il[i] = 0; 105978910aadSHong Zhang } 106078910aadSHong Zhang 106178910aadSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 106278910aadSHong Zhang nlnk = am + 1; 106378910aadSHong Zhang ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 106478910aadSHong Zhang 106595b5ac22SHong Zhang /* initial FreeSpace size is fill*(ai[am]+am)/2 */ 106695b5ac22SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+am)/2),&free_space);CHKERRQ(ierr); 106726fbe8dcSKarl Rupp 106878910aadSHong Zhang current_space = free_space; 106926fbe8dcSKarl Rupp 107095b5ac22SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+am)/2),&free_space_lvl);CHKERRQ(ierr); 107178910aadSHong Zhang current_space_lvl = free_space_lvl; 107278910aadSHong Zhang 107378910aadSHong Zhang for (k=0; k<am; k++) { /* for each active row k */ 107478910aadSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 107578910aadSHong Zhang nzk = 0; 107678910aadSHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 107778910aadSHong Zhang ncols_upper = 0; 107878910aadSHong Zhang cols = cols_lvl + am; 107978910aadSHong Zhang for (j=0; j<ncols; j++) { 108078910aadSHong Zhang i = rip[*(aj + ai[rip[k]] + j)]; 108178910aadSHong Zhang if (i >= k) { /* only take upper triangular entry */ 108278910aadSHong Zhang cols[ncols_upper] = i; 108378910aadSHong Zhang cols_lvl[ncols_upper] = -1; /* initialize level for nonzero entries */ 108478910aadSHong Zhang ncols_upper++; 108578910aadSHong Zhang } 108678910aadSHong Zhang } 108778910aadSHong Zhang ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 108878910aadSHong Zhang nzk += nlnk; 108978910aadSHong Zhang 109078910aadSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 109178910aadSHong Zhang prow = jl[k]; /* 1st pivot row */ 109278910aadSHong Zhang 109378910aadSHong Zhang while (prow < k) { 109478910aadSHong Zhang nextprow = jl[prow]; 109578910aadSHong Zhang 109678910aadSHong Zhang /* merge prow into k-th row */ 109778910aadSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 109878910aadSHong Zhang jmax = ui[prow+1]; 109978910aadSHong Zhang ncols = jmax-jmin; 110078910aadSHong Zhang i = jmin - ui[prow]; 110178910aadSHong Zhang cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 110278910aadSHong Zhang for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j); 11035a8e39fbSHong Zhang ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 110478910aadSHong Zhang nzk += nlnk; 110578910aadSHong Zhang 110678910aadSHong Zhang /* update il and jl for prow */ 110778910aadSHong Zhang if (jmin < jmax) { 110878910aadSHong Zhang il[prow] = jmin; 110926fbe8dcSKarl Rupp 111078910aadSHong Zhang j = *cols; jl[prow] = jl[j]; jl[j] = prow; 111178910aadSHong Zhang } 111278910aadSHong Zhang prow = nextprow; 111378910aadSHong Zhang } 111478910aadSHong Zhang 111578910aadSHong Zhang /* if free space is not available, make more free space */ 111678910aadSHong Zhang if (current_space->local_remaining<nzk) { 111778910aadSHong Zhang i = am - k + 1; /* num of unfactored rows */ 111878910aadSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 1119a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 1120a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space_lvl);CHKERRQ(ierr); 112178910aadSHong Zhang reallocs++; 112278910aadSHong Zhang } 112378910aadSHong Zhang 112478910aadSHong Zhang /* copy data into free_space and free_space_lvl, then initialize lnk */ 112578910aadSHong Zhang ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 112678910aadSHong Zhang 112778910aadSHong Zhang /* add the k-th row into il and jl */ 112878910aadSHong Zhang if (nzk-1 > 0) { 112978910aadSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 113078910aadSHong Zhang jl[k] = jl[i]; jl[i] = k; 113178910aadSHong Zhang il[k] = ui[k] + 1; 113278910aadSHong Zhang } 113378910aadSHong Zhang uj_ptr[k] = current_space->array; 113478910aadSHong Zhang uj_lvl_ptr[k] = current_space_lvl->array; 113578910aadSHong Zhang 113678910aadSHong Zhang current_space->array += nzk; 113778910aadSHong Zhang current_space->local_used += nzk; 113878910aadSHong Zhang current_space->local_remaining -= nzk; 113978910aadSHong Zhang 114078910aadSHong Zhang current_space_lvl->array += nzk; 114178910aadSHong Zhang current_space_lvl->local_used += nzk; 114278910aadSHong Zhang current_space_lvl->local_remaining -= nzk; 114378910aadSHong Zhang 114478910aadSHong Zhang ui[k+1] = ui[k] + nzk; 114578910aadSHong Zhang } 114678910aadSHong Zhang 114778910aadSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 1148fca92195SBarry Smith ierr = PetscFree4(uj_ptr,uj_lvl_ptr,il,jl);CHKERRQ(ierr); 114978910aadSHong Zhang ierr = PetscFree(cols_lvl);CHKERRQ(ierr); 115078910aadSHong Zhang 11519263d837SHong Zhang /* copy free_space into uj and free free_space; set uj in new datastructure; */ 1152785e854fSJed Brown ierr = PetscMalloc1((ui[am]+1),&uj);CHKERRQ(ierr); 1153a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 115478910aadSHong Zhang ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1155a1a86e44SBarry Smith ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 115678910aadSHong Zhang 115778910aadSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 1158719d5645SBarry Smith B = fact; 11590298fd71SBarry Smith ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); 116078910aadSHong Zhang 116178910aadSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 116278910aadSHong Zhang b->singlemalloc = PETSC_FALSE; 1163e6b907acSBarry Smith b->free_a = PETSC_TRUE; 1164e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 116526fbe8dcSKarl Rupp 1166785e854fSJed Brown ierr = PetscMalloc1((ui[am]+1),&b->a);CHKERRQ(ierr); 116726fbe8dcSKarl Rupp 116878910aadSHong Zhang b->j = uj; 116978910aadSHong Zhang b->i = ui; 117078910aadSHong Zhang b->diag = 0; 117178910aadSHong Zhang b->ilen = 0; 117278910aadSHong Zhang b->imax = 0; 117378910aadSHong Zhang b->row = perm; 117478910aadSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 117526fbe8dcSKarl Rupp 117678910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 117726fbe8dcSKarl Rupp 117878910aadSHong Zhang b->icol = perm; 117926fbe8dcSKarl Rupp 118078910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1181785e854fSJed Brown ierr = PetscMalloc1((am+1),&b->solve_work);CHKERRQ(ierr); 11823bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 118326fbe8dcSKarl Rupp 118478910aadSHong Zhang b->maxnz = b->nz = ui[am]; 118578910aadSHong Zhang 118678910aadSHong Zhang B->info.factor_mallocs = reallocs; 118778910aadSHong Zhang B->info.fill_ratio_given = fill; 118875567043SBarry Smith if (ai[am] != 0.) { 118995b5ac22SHong Zhang /* nonzeros in lower triangular part of A (includign diagonals)= (ai[am]+am)/2 */ 119095b5ac22SHong Zhang B->info.fill_ratio_needed = ((PetscReal)2*ui[am])/(ai[am]+am); 119178910aadSHong Zhang } else { 119278910aadSHong Zhang B->info.fill_ratio_needed = 0.0; 119378910aadSHong Zhang } 11949263d837SHong Zhang #if defined(PETSC_USE_INFO) 11959263d837SHong Zhang if (ai[am] != 0) { 119695b5ac22SHong Zhang PetscReal af = B->info.fill_ratio_needed; 119757622a8eSBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);CHKERRQ(ierr); 119857622a8eSBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr); 119957622a8eSBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);CHKERRQ(ierr); 12009263d837SHong Zhang } else { 12019263d837SHong Zhang ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 12029263d837SHong Zhang } 12039263d837SHong Zhang #endif 120478910aadSHong Zhang if (perm_identity) { 12050a3c351aSHong Zhang B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 12060a3c351aSHong Zhang B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 120778910aadSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 120878910aadSHong Zhang } else { 1209719d5645SBarry Smith (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 121078910aadSHong Zhang } 1211c05c3958SHong Zhang PetscFunctionReturn(0); 1212c05c3958SHong Zhang } 1213c05c3958SHong Zhang 1214c05c3958SHong Zhang #undef __FUNCT__ 1215c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ" 12160481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 1217c05c3958SHong Zhang { 121878910aadSHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 121978910aadSHong Zhang Mat_SeqSBAIJ *b; 122078910aadSHong Zhang Mat B; 122178910aadSHong Zhang PetscErrorCode ierr; 12229186b105SHong Zhang PetscBool perm_identity,missing; 122378910aadSHong Zhang PetscReal fill = info->fill; 12245d0c19d7SBarry Smith const PetscInt *rip; 12255d0c19d7SBarry Smith PetscInt i,mbs=a->mbs,bs=A->rmap->bs,*ai=a->i,*aj=a->j,reallocs=0,prow; 122678910aadSHong Zhang PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 122778910aadSHong Zhang PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 12280298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 122978910aadSHong Zhang PetscBT lnkbt; 123078910aadSHong Zhang 1231c05c3958SHong Zhang PetscFunctionBegin; 12326ad2eaddSHong Zhang if (bs > 1) { /* convert to seqsbaij */ 12336ad2eaddSHong Zhang if (!a->sbaijMat) { 1234ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 12356ad2eaddSHong Zhang } 1236719d5645SBarry Smith (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */ 123726fbe8dcSKarl Rupp 1238719d5645SBarry Smith ierr = MatCholeskyFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr); 12396ad2eaddSHong Zhang PetscFunctionReturn(0); 12406ad2eaddSHong Zhang } 12416ad2eaddSHong Zhang 12429186b105SHong Zhang ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr); 12439186b105SHong Zhang if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i); 12449186b105SHong Zhang 124578910aadSHong Zhang /* check whether perm is the identity mapping */ 124678910aadSHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 1247e32f2f54SBarry Smith if (!perm_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported"); 124878910aadSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 124978910aadSHong Zhang 125078910aadSHong Zhang /* initialization */ 1251785e854fSJed Brown ierr = PetscMalloc1((mbs+1),&ui);CHKERRQ(ierr); 1252e60cf9a0SBarry Smith ui[0] = 0; 125378910aadSHong Zhang 125478910aadSHong Zhang /* jl: linked list for storing indices of the pivot rows 125578910aadSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */ 1256dcca6d9dSJed Brown ierr = PetscMalloc4(mbs,&ui_ptr,mbs,&il,mbs,&jl,mbs,&cols);CHKERRQ(ierr); 125778910aadSHong Zhang for (i=0; i<mbs; i++) { 1258e60cf9a0SBarry Smith jl[i] = mbs; il[i] = 0; 125978910aadSHong Zhang } 126078910aadSHong Zhang 126178910aadSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 126278910aadSHong Zhang nlnk = mbs + 1; 126378910aadSHong Zhang ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 126478910aadSHong Zhang 12656a69fef8SHong Zhang /* initial FreeSpace size is fill* (ai[mbs]+mbs)/2 */ 126695b5ac22SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+mbs)/2),&free_space);CHKERRQ(ierr); 126726fbe8dcSKarl Rupp 126878910aadSHong Zhang current_space = free_space; 126978910aadSHong Zhang 127078910aadSHong Zhang for (k=0; k<mbs; k++) { /* for each active row k */ 127178910aadSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 127278910aadSHong Zhang nzk = 0; 127378910aadSHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 1274e60cf9a0SBarry Smith ncols_upper = 0; 127578910aadSHong Zhang for (j=0; j<ncols; j++) { 127678910aadSHong Zhang i = rip[*(aj + ai[rip[k]] + j)]; 127778910aadSHong Zhang if (i >= k) { /* only take upper triangular entry */ 127878910aadSHong Zhang cols[ncols_upper] = i; 127978910aadSHong Zhang ncols_upper++; 128078910aadSHong Zhang } 128178910aadSHong Zhang } 128278910aadSHong Zhang ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 128378910aadSHong Zhang nzk += nlnk; 128478910aadSHong Zhang 128578910aadSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 128678910aadSHong Zhang prow = jl[k]; /* 1st pivot row */ 128778910aadSHong Zhang 128878910aadSHong Zhang while (prow < k) { 128978910aadSHong Zhang nextprow = jl[prow]; 129078910aadSHong Zhang /* merge prow into k-th row */ 129178910aadSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */ 129278910aadSHong Zhang jmax = ui[prow+1]; 129378910aadSHong Zhang ncols = jmax-jmin; 129478910aadSHong Zhang uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */ 12955a8e39fbSHong Zhang ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 129678910aadSHong Zhang nzk += nlnk; 129778910aadSHong Zhang 129878910aadSHong Zhang /* update il and jl for prow */ 129978910aadSHong Zhang if (jmin < jmax) { 130078910aadSHong Zhang il[prow] = jmin; 130126fbe8dcSKarl Rupp j = *uj_ptr; 130226fbe8dcSKarl Rupp jl[prow] = jl[j]; 130326fbe8dcSKarl Rupp jl[j] = prow; 130478910aadSHong Zhang } 130578910aadSHong Zhang prow = nextprow; 130678910aadSHong Zhang } 130778910aadSHong Zhang 130878910aadSHong Zhang /* if free space is not available, make more free space */ 130978910aadSHong Zhang if (current_space->local_remaining<nzk) { 131078910aadSHong Zhang i = mbs - k + 1; /* num of unfactored rows */ 131178910aadSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 1312a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 131378910aadSHong Zhang reallocs++; 131478910aadSHong Zhang } 131578910aadSHong Zhang 131678910aadSHong Zhang /* copy data into free space, then initialize lnk */ 131778910aadSHong Zhang ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 131878910aadSHong Zhang 131978910aadSHong Zhang /* add the k-th row into il and jl */ 132078910aadSHong Zhang if (nzk-1 > 0) { 132178910aadSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ 132278910aadSHong Zhang jl[k] = jl[i]; jl[i] = k; 132378910aadSHong Zhang il[k] = ui[k] + 1; 132478910aadSHong Zhang } 132578910aadSHong Zhang ui_ptr[k] = current_space->array; 132678910aadSHong Zhang current_space->array += nzk; 132778910aadSHong Zhang current_space->local_used += nzk; 132878910aadSHong Zhang current_space->local_remaining -= nzk; 132978910aadSHong Zhang 133078910aadSHong Zhang ui[k+1] = ui[k] + nzk; 133178910aadSHong Zhang } 133278910aadSHong Zhang 133378910aadSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 1334fca92195SBarry Smith ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr); 133578910aadSHong Zhang 13369263d837SHong Zhang /* copy free_space into uj and free free_space; set uj in new datastructure; */ 1337785e854fSJed Brown ierr = PetscMalloc1((ui[mbs]+1),&uj);CHKERRQ(ierr); 1338a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 133978910aadSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 134078910aadSHong Zhang 134178910aadSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 1342719d5645SBarry Smith B = fact; 13430298fd71SBarry Smith ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); 134478910aadSHong Zhang 134578910aadSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 134678910aadSHong Zhang b->singlemalloc = PETSC_FALSE; 1347e6b907acSBarry Smith b->free_a = PETSC_TRUE; 1348e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 134926fbe8dcSKarl Rupp 1350785e854fSJed Brown ierr = PetscMalloc1((ui[mbs]+1),&b->a);CHKERRQ(ierr); 135126fbe8dcSKarl Rupp 135278910aadSHong Zhang b->j = uj; 135378910aadSHong Zhang b->i = ui; 135478910aadSHong Zhang b->diag = 0; 135578910aadSHong Zhang b->ilen = 0; 135678910aadSHong Zhang b->imax = 0; 135778910aadSHong Zhang b->row = perm; 135878910aadSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 135926fbe8dcSKarl Rupp 136078910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 136178910aadSHong Zhang b->icol = perm; 136278910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1363785e854fSJed Brown ierr = PetscMalloc1((mbs+1),&b->solve_work);CHKERRQ(ierr); 13643bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 136578910aadSHong Zhang b->maxnz = b->nz = ui[mbs]; 136678910aadSHong Zhang 136778910aadSHong Zhang B->info.factor_mallocs = reallocs; 136878910aadSHong Zhang B->info.fill_ratio_given = fill; 136975567043SBarry Smith if (ai[mbs] != 0.) { 137095b5ac22SHong Zhang /* nonzeros in lower triangular part of A = (ai[mbs]+mbs)/2 */ 137195b5ac22SHong Zhang B->info.fill_ratio_needed = ((PetscReal)2*ui[mbs])/(ai[mbs]+mbs); 137278910aadSHong Zhang } else { 137378910aadSHong Zhang B->info.fill_ratio_needed = 0.0; 137478910aadSHong Zhang } 13759263d837SHong Zhang #if defined(PETSC_USE_INFO) 13769263d837SHong Zhang if (ai[mbs] != 0.) { 13779263d837SHong Zhang PetscReal af = B->info.fill_ratio_needed; 137857622a8eSBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);CHKERRQ(ierr); 137957622a8eSBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr); 138057622a8eSBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);CHKERRQ(ierr); 13819263d837SHong Zhang } else { 13829263d837SHong Zhang ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 13839263d837SHong Zhang } 13849263d837SHong Zhang #endif 138578910aadSHong Zhang if (perm_identity) { 13866ad2eaddSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 138778910aadSHong Zhang } else { 13886ad2eaddSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 138978910aadSHong Zhang } 1390c05c3958SHong Zhang PetscFunctionReturn(0); 1391c05c3958SHong Zhang } 1392c8342467SHong Zhang 13932efa7f71SHong Zhang #undef __FUNCT__ 13944dd39f65SShri Abhyankar #define __FUNCT__ "MatSolve_SeqBAIJ_N_NaturalOrdering" 13954dd39f65SShri Abhyankar PetscErrorCode MatSolve_SeqBAIJ_N_NaturalOrdering(Mat A,Vec bb,Vec xx) 13961a83e813SShri Abhyankar { 13971a83e813SShri Abhyankar Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 13981a83e813SShri Abhyankar PetscErrorCode ierr; 13991a83e813SShri Abhyankar const PetscInt *ai=a->i,*aj=a->j,*adiag=a->diag,*vi; 14001a83e813SShri Abhyankar PetscInt i,k,n=a->mbs; 14011a83e813SShri Abhyankar PetscInt nz,bs=A->rmap->bs,bs2=a->bs2; 14021a83e813SShri Abhyankar MatScalar *aa=a->a,*v; 14031a83e813SShri Abhyankar PetscScalar *x,*b,*s,*t,*ls; 14041a83e813SShri Abhyankar 14051a83e813SShri Abhyankar PetscFunctionBegin; 14061a83e813SShri Abhyankar ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 14071a83e813SShri Abhyankar ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 14081a83e813SShri Abhyankar t = a->solve_work; 14091a83e813SShri Abhyankar 14101a83e813SShri Abhyankar /* forward solve the lower triangular */ 14111a83e813SShri Abhyankar ierr = PetscMemcpy(t,b,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy 1st block of b to t */ 14121a83e813SShri Abhyankar 14131a83e813SShri Abhyankar for (i=1; i<n; i++) { 14141a83e813SShri Abhyankar v = aa + bs2*ai[i]; 14151a83e813SShri Abhyankar vi = aj + ai[i]; 14161a83e813SShri Abhyankar nz = ai[i+1] - ai[i]; 14171a83e813SShri Abhyankar s = t + bs*i; 14181a83e813SShri Abhyankar ierr = PetscMemcpy(s,b+bs*i,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy i_th block of b to t */ 14191a83e813SShri Abhyankar for (k=0;k<nz;k++) { 142096b95a6bSBarry Smith PetscKernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[k]); 14211a83e813SShri Abhyankar v += bs2; 14221a83e813SShri Abhyankar } 14231a83e813SShri Abhyankar } 14241a83e813SShri Abhyankar 14251a83e813SShri Abhyankar /* backward solve the upper triangular */ 14261a83e813SShri Abhyankar ls = a->solve_work + A->cmap->n; 14271a83e813SShri Abhyankar for (i=n-1; i>=0; i--) { 14281a83e813SShri Abhyankar v = aa + bs2*(adiag[i+1]+1); 14291a83e813SShri Abhyankar vi = aj + adiag[i+1]+1; 14301a83e813SShri Abhyankar nz = adiag[i] - adiag[i+1]-1; 14311a83e813SShri Abhyankar ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 14321a83e813SShri Abhyankar for (k=0; k<nz; k++) { 143396b95a6bSBarry Smith PetscKernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[k]); 14341a83e813SShri Abhyankar v += bs2; 14351a83e813SShri Abhyankar } 143696b95a6bSBarry Smith PetscKernel_w_gets_A_times_v(bs,ls,aa+bs2*adiag[i],t+i*bs); /* *inv(diagonal[i]) */ 14371a83e813SShri Abhyankar ierr = PetscMemcpy(x+i*bs,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 14381a83e813SShri Abhyankar } 14391a83e813SShri Abhyankar 14401a83e813SShri Abhyankar ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 14411a83e813SShri Abhyankar ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 14421a83e813SShri Abhyankar ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr); 14431a83e813SShri Abhyankar PetscFunctionReturn(0); 14441a83e813SShri Abhyankar } 14451a83e813SShri Abhyankar 14462efa7f71SHong Zhang #undef __FUNCT__ 14474dd39f65SShri Abhyankar #define __FUNCT__ "MatSolve_SeqBAIJ_N" 14484dd39f65SShri Abhyankar PetscErrorCode MatSolve_SeqBAIJ_N(Mat A,Vec bb,Vec xx) 144935aa4fcfSShri Abhyankar { 145035aa4fcfSShri Abhyankar Mat_SeqBAIJ *a =(Mat_SeqBAIJ*)A->data; 145135aa4fcfSShri Abhyankar IS iscol=a->col,isrow=a->row; 145235aa4fcfSShri Abhyankar PetscErrorCode ierr; 145335aa4fcfSShri Abhyankar const PetscInt *r,*c,*rout,*cout,*ai=a->i,*aj=a->j,*adiag=a->diag,*vi; 145435aa4fcfSShri Abhyankar PetscInt i,m,n=a->mbs; 145535aa4fcfSShri Abhyankar PetscInt nz,bs=A->rmap->bs,bs2=a->bs2; 145635aa4fcfSShri Abhyankar MatScalar *aa=a->a,*v; 145735aa4fcfSShri Abhyankar PetscScalar *x,*b,*s,*t,*ls; 145835aa4fcfSShri Abhyankar 145935aa4fcfSShri Abhyankar PetscFunctionBegin; 146035aa4fcfSShri Abhyankar ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 146135aa4fcfSShri Abhyankar ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 146235aa4fcfSShri Abhyankar t = a->solve_work; 146335aa4fcfSShri Abhyankar 146435aa4fcfSShri Abhyankar ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 146535aa4fcfSShri Abhyankar ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 146635aa4fcfSShri Abhyankar 146735aa4fcfSShri Abhyankar /* forward solve the lower triangular */ 146835aa4fcfSShri Abhyankar ierr = PetscMemcpy(t,b+bs*r[0],bs*sizeof(PetscScalar));CHKERRQ(ierr); 146935aa4fcfSShri Abhyankar for (i=1; i<n; i++) { 147035aa4fcfSShri Abhyankar v = aa + bs2*ai[i]; 147135aa4fcfSShri Abhyankar vi = aj + ai[i]; 147235aa4fcfSShri Abhyankar nz = ai[i+1] - ai[i]; 147335aa4fcfSShri Abhyankar s = t + bs*i; 147435aa4fcfSShri Abhyankar ierr = PetscMemcpy(s,b+bs*r[i],bs*sizeof(PetscScalar));CHKERRQ(ierr); 147535aa4fcfSShri Abhyankar for (m=0; m<nz; m++) { 147696b95a6bSBarry Smith PetscKernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[m]); 147735aa4fcfSShri Abhyankar v += bs2; 147835aa4fcfSShri Abhyankar } 147935aa4fcfSShri Abhyankar } 148035aa4fcfSShri Abhyankar 148135aa4fcfSShri Abhyankar /* backward solve the upper triangular */ 148235aa4fcfSShri Abhyankar ls = a->solve_work + A->cmap->n; 148335aa4fcfSShri Abhyankar for (i=n-1; i>=0; i--) { 148435aa4fcfSShri Abhyankar v = aa + bs2*(adiag[i+1]+1); 148535aa4fcfSShri Abhyankar vi = aj + adiag[i+1]+1; 148635aa4fcfSShri Abhyankar nz = adiag[i] - adiag[i+1] - 1; 148735aa4fcfSShri Abhyankar ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 148835aa4fcfSShri Abhyankar for (m=0; m<nz; m++) { 148996b95a6bSBarry Smith PetscKernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[m]); 149035aa4fcfSShri Abhyankar v += bs2; 149135aa4fcfSShri Abhyankar } 149296b95a6bSBarry Smith PetscKernel_w_gets_A_times_v(bs,ls,v,t+i*bs); /* *inv(diagonal[i]) */ 149335aa4fcfSShri Abhyankar ierr = PetscMemcpy(x + bs*c[i],t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 149435aa4fcfSShri Abhyankar } 149535aa4fcfSShri Abhyankar ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 149635aa4fcfSShri Abhyankar ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 149735aa4fcfSShri Abhyankar ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 149835aa4fcfSShri Abhyankar ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 149935aa4fcfSShri Abhyankar ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr); 150035aa4fcfSShri Abhyankar PetscFunctionReturn(0); 150135aa4fcfSShri Abhyankar } 150235aa4fcfSShri Abhyankar 150335aa4fcfSShri Abhyankar #undef __FUNCT__ 1504a25532f0SBarry Smith #define __FUNCT__ "MatBlockAbs_privat" 1505a25532f0SBarry Smith /* 1506a25532f0SBarry Smith For each block in an block array saves the largest absolute value in the block into another array 1507a25532f0SBarry Smith */ 1508a25532f0SBarry Smith static PetscErrorCode MatBlockAbs_private(PetscInt nbs,PetscInt bs2,PetscScalar *blockarray,PetscReal *absarray) 15092efa7f71SHong Zhang { 15102efa7f71SHong Zhang PetscErrorCode ierr; 15112efa7f71SHong Zhang PetscInt i,j; 15125fd66863SKarl Rupp 15132efa7f71SHong Zhang PetscFunctionBegin; 15142efa7f71SHong Zhang ierr = PetscMemzero(absarray,(nbs+1)*sizeof(PetscReal));CHKERRQ(ierr); 15152efa7f71SHong Zhang for (i=0; i<nbs; i++) { 15162efa7f71SHong Zhang for (j=0; j<bs2; j++) { 15172efa7f71SHong Zhang if (absarray[i] < PetscAbsScalar(blockarray[i*nbs+j])) absarray[i] = PetscAbsScalar(blockarray[i*nbs+j]); 15182efa7f71SHong Zhang } 15192efa7f71SHong Zhang } 15202efa7f71SHong Zhang PetscFunctionReturn(0); 15212efa7f71SHong Zhang } 15222efa7f71SHong Zhang 1523c8342467SHong Zhang #undef __FUNCT__ 1524c8342467SHong Zhang #define __FUNCT__ "MatILUDTFactor_SeqBAIJ" 1525fe97e370SBarry Smith /* 1526fe97e370SBarry Smith This needs to be renamed and called by the regular MatILUFactor_SeqBAIJ when drop tolerance is used 1527fe97e370SBarry Smith */ 1528c8342467SHong Zhang PetscErrorCode MatILUDTFactor_SeqBAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact) 1529c8342467SHong Zhang { 15302efa7f71SHong Zhang Mat B = *fact; 15312efa7f71SHong Zhang Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b; 15322efa7f71SHong Zhang IS isicol; 15332efa7f71SHong Zhang PetscErrorCode ierr; 15342efa7f71SHong Zhang const PetscInt *r,*ic; 15352efa7f71SHong Zhang PetscInt i,mbs=a->mbs,bs=A->rmap->bs,bs2=a->bs2,*ai=a->i,*aj=a->j,*ajtmp,*adiag; 15362efa7f71SHong Zhang PetscInt *bi,*bj,*bdiag; 15372efa7f71SHong Zhang 15382efa7f71SHong Zhang PetscInt row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au; 15392efa7f71SHong Zhang PetscInt nlnk,*lnk; 15402efa7f71SHong Zhang PetscBT lnkbt; 1541ace3abfcSBarry Smith PetscBool row_identity,icol_identity; 15422efa7f71SHong Zhang MatScalar *aatmp,*pv,*batmp,*ba,*rtmp,*pc,*multiplier,*vtmp; 15432efa7f71SHong Zhang PetscInt j,nz,*pj,*bjtmp,k,ncut,*jtmp; 15442efa7f71SHong Zhang 1545182b8fbaSHong Zhang PetscReal dt=info->dt; /* shift=info->shiftamount; */ 15462efa7f71SHong Zhang PetscInt nnz_max; 1547ace3abfcSBarry Smith PetscBool missing; 1548021822bcSHong Zhang PetscReal *vtmp_abs; 154997ef94ebSSatish Balay MatScalar *v_work; 155097ef94ebSSatish Balay PetscInt *v_pivots; 15512efa7f71SHong Zhang 1552c8342467SHong Zhang PetscFunctionBegin; 15532efa7f71SHong Zhang /* ------- symbolic factorization, can be reused ---------*/ 15542efa7f71SHong Zhang ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr); 1555e32f2f54SBarry Smith if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i); 15562efa7f71SHong Zhang adiag=a->diag; 15572efa7f71SHong Zhang 15582efa7f71SHong Zhang ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 15592efa7f71SHong Zhang 15602efa7f71SHong Zhang /* bdiag is location of diagonal in factor */ 1561785e854fSJed Brown ierr = PetscMalloc1((mbs+1),&bdiag);CHKERRQ(ierr); 15622efa7f71SHong Zhang 15632efa7f71SHong Zhang /* allocate row pointers bi */ 1564785e854fSJed Brown ierr = PetscMalloc1((2*mbs+2),&bi);CHKERRQ(ierr); 15652efa7f71SHong Zhang 15662efa7f71SHong Zhang /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */ 15672efa7f71SHong Zhang dtcount = (PetscInt)info->dtcount; 15686bce7ff8SHong Zhang if (dtcount > mbs-1) dtcount = mbs-1; 15696bce7ff8SHong Zhang nnz_max = ai[mbs]+2*mbs*dtcount +2; 15706da515a1SHong Zhang /* printf("MatILUDTFactor_SeqBAIJ, bs %d, ai[mbs] %d, nnz_max %d, dtcount %d\n",bs,ai[mbs],nnz_max,dtcount); */ 1571785e854fSJed Brown ierr = PetscMalloc1(nnz_max,&bj);CHKERRQ(ierr); 15726bce7ff8SHong Zhang nnz_max = nnz_max*bs2; 1573785e854fSJed Brown ierr = PetscMalloc1(nnz_max,&ba);CHKERRQ(ierr); 15742efa7f71SHong Zhang 15752efa7f71SHong Zhang /* put together the new matrix */ 15760298fd71SBarry Smith ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); 15773bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)B,(PetscObject)isicol);CHKERRQ(ierr); 157826fbe8dcSKarl Rupp 15792efa7f71SHong Zhang b = (Mat_SeqBAIJ*)(B)->data; 15802efa7f71SHong Zhang b->free_a = PETSC_TRUE; 15812efa7f71SHong Zhang b->free_ij = PETSC_TRUE; 15822efa7f71SHong Zhang b->singlemalloc = PETSC_FALSE; 158326fbe8dcSKarl Rupp 15842efa7f71SHong Zhang b->a = ba; 15852efa7f71SHong Zhang b->j = bj; 15862efa7f71SHong Zhang b->i = bi; 15872efa7f71SHong Zhang b->diag = bdiag; 15882efa7f71SHong Zhang b->ilen = 0; 15892efa7f71SHong Zhang b->imax = 0; 15902efa7f71SHong Zhang b->row = isrow; 15912efa7f71SHong Zhang b->col = iscol; 159226fbe8dcSKarl Rupp 15932efa7f71SHong Zhang ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 15942efa7f71SHong Zhang ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 159526fbe8dcSKarl Rupp 15962efa7f71SHong Zhang b->icol = isicol; 1597785e854fSJed Brown ierr = PetscMalloc1((bs*(mbs+1)),&b->solve_work);CHKERRQ(ierr); 15983bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 15992efa7f71SHong Zhang b->maxnz = nnz_max/bs2; 16002efa7f71SHong Zhang 1601d5f3da31SBarry Smith (B)->factortype = MAT_FACTOR_ILUDT; 16022efa7f71SHong Zhang (B)->info.factor_mallocs = 0; 16032efa7f71SHong Zhang (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)(ai[mbs]*bs2)); 16042efa7f71SHong Zhang /* ------- end of symbolic factorization ---------*/ 16052efa7f71SHong Zhang ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 16062efa7f71SHong Zhang ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 16072efa7f71SHong Zhang 16082efa7f71SHong Zhang /* linked list for storing column indices of the active row */ 16092efa7f71SHong Zhang nlnk = mbs + 1; 16102efa7f71SHong Zhang ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 16112efa7f71SHong Zhang 16122efa7f71SHong Zhang /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */ 1613dcca6d9dSJed Brown ierr = PetscMalloc2(mbs,&im,mbs,&jtmp);CHKERRQ(ierr); 16142efa7f71SHong Zhang /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */ 1615dcca6d9dSJed Brown ierr = PetscMalloc2(mbs*bs2,&rtmp,mbs*bs2,&vtmp);CHKERRQ(ierr); 1616785e854fSJed Brown ierr = PetscMalloc1((mbs+1),&vtmp_abs);CHKERRQ(ierr); 1617dcca6d9dSJed Brown ierr = PetscMalloc3(bs,&v_work,bs2,&multiplier,bs,&v_pivots);CHKERRQ(ierr); 16182efa7f71SHong Zhang 1619e60cf9a0SBarry Smith bi[0] = 0; 16202efa7f71SHong Zhang bdiag[0] = (nnz_max/bs2)-1; /* location of diagonal in factor B */ 16216bce7ff8SHong Zhang bi[2*mbs+1] = bdiag[0]+1; /* endof bj and ba array */ 16222efa7f71SHong Zhang for (i=0; i<mbs; i++) { 16232efa7f71SHong Zhang /* copy initial fill into linked list */ 16242efa7f71SHong Zhang nzi = 0; /* nonzeros for active row i */ 16252efa7f71SHong Zhang nzi = ai[r[i]+1] - ai[r[i]]; 1626e32f2f54SBarry 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); 16272efa7f71SHong Zhang nzi_al = adiag[r[i]] - ai[r[i]]; 16282efa7f71SHong Zhang nzi_au = ai[r[i]+1] - adiag[r[i]] -1; 16296da515a1SHong Zhang /* printf("row %d, nzi_al/au %d %d\n",i,nzi_al,nzi_au); */ 16302efa7f71SHong Zhang 16312efa7f71SHong Zhang /* load in initial unfactored row */ 16322efa7f71SHong Zhang ajtmp = aj + ai[r[i]]; 16332efa7f71SHong Zhang ierr = PetscLLAddPerm(nzi,ajtmp,ic,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1634fca92195SBarry Smith ierr = PetscMemzero(rtmp,mbs*bs2*sizeof(PetscScalar));CHKERRQ(ierr); 16352efa7f71SHong Zhang aatmp = a->a + bs2*ai[r[i]]; 16362efa7f71SHong Zhang for (j=0; j<nzi; j++) { 16372efa7f71SHong Zhang ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],aatmp+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 16382efa7f71SHong Zhang } 16392efa7f71SHong Zhang 16402efa7f71SHong Zhang /* add pivot rows into linked list */ 16412efa7f71SHong Zhang row = lnk[mbs]; 16422efa7f71SHong Zhang while (row < i) { 16432efa7f71SHong Zhang nzi_bl = bi[row+1] - bi[row] + 1; 16442efa7f71SHong Zhang bjtmp = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */ 16452efa7f71SHong Zhang ierr = PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);CHKERRQ(ierr); 16462efa7f71SHong Zhang nzi += nlnk; 16472efa7f71SHong Zhang row = lnk[row]; 16482efa7f71SHong Zhang } 16492efa7f71SHong Zhang 16502efa7f71SHong Zhang /* copy data from lnk into jtmp, then initialize lnk */ 16512efa7f71SHong Zhang ierr = PetscLLClean(mbs,mbs,nzi,lnk,jtmp,lnkbt);CHKERRQ(ierr); 16522efa7f71SHong Zhang 16532efa7f71SHong Zhang /* numerical factorization */ 16542efa7f71SHong Zhang bjtmp = jtmp; 16552efa7f71SHong Zhang row = *bjtmp++; /* 1st pivot row */ 16562efa7f71SHong Zhang 16572efa7f71SHong Zhang while (row < i) { 16582efa7f71SHong Zhang pc = rtmp + bs2*row; 16592efa7f71SHong Zhang pv = ba + bs2*bdiag[row]; /* inv(diag) of the pivot row */ 166096b95a6bSBarry Smith PetscKernel_A_gets_A_times_B(bs,pc,pv,multiplier); /* pc= multiplier = pc*inv(diag[row]) */ 1661a25532f0SBarry Smith ierr = MatBlockAbs_private(1,bs2,pc,vtmp_abs);CHKERRQ(ierr); 16622efa7f71SHong Zhang if (vtmp_abs[0] > dt) { /* apply tolerance dropping rule */ 16632efa7f71SHong Zhang pj = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */ 16642efa7f71SHong Zhang pv = ba + bs2*(bdiag[row+1] + 1); 16652efa7f71SHong Zhang nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */ 16662efa7f71SHong Zhang for (j=0; j<nz; j++) { 166796b95a6bSBarry Smith PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); 16682efa7f71SHong Zhang } 16696da515a1SHong Zhang /* ierr = PetscLogFlops(bslog*(nz+1.0)-bs);CHKERRQ(ierr); */ 16702efa7f71SHong Zhang } 16712efa7f71SHong Zhang row = *bjtmp++; 16722efa7f71SHong Zhang } 16732efa7f71SHong Zhang 16742efa7f71SHong Zhang /* copy sparse rtmp into contiguous vtmp; separate L and U part */ 16752efa7f71SHong Zhang nzi_bl = 0; j = 0; 16762efa7f71SHong Zhang while (jtmp[j] < i) { /* L-part. Note: jtmp is sorted */ 16772efa7f71SHong Zhang ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 16782efa7f71SHong Zhang nzi_bl++; j++; 16792efa7f71SHong Zhang } 16802efa7f71SHong Zhang nzi_bu = nzi - nzi_bl -1; 16816da515a1SHong Zhang /* printf("nzi %d, nzi_bl %d, nzi_bu %d\n",nzi,nzi_bl,nzi_bu); */ 16822efa7f71SHong Zhang 16832efa7f71SHong Zhang while (j < nzi) { /* U-part */ 16842efa7f71SHong Zhang ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 16852efa7f71SHong Zhang /* 16862efa7f71SHong Zhang printf(" col %d: ",jtmp[j]); 16872efa7f71SHong Zhang for (j1=0; j1<bs2; j1++) printf(" %g",*(vtmp+bs2*j+j1)); 16882efa7f71SHong Zhang printf(" \n"); 16892efa7f71SHong Zhang */ 16902efa7f71SHong Zhang j++; 16912efa7f71SHong Zhang } 16922efa7f71SHong Zhang 1693a25532f0SBarry Smith ierr = MatBlockAbs_private(nzi,bs2,vtmp,vtmp_abs);CHKERRQ(ierr); 16942efa7f71SHong Zhang /* 16952efa7f71SHong Zhang printf(" row %d, nzi %d, vtmp_abs\n",i,nzi); 16962efa7f71SHong Zhang for (j1=0; j1<nzi; j1++) printf(" (%d %g),",jtmp[j1],vtmp_abs[j1]); 16972efa7f71SHong Zhang printf(" \n"); 16982efa7f71SHong Zhang */ 16992efa7f71SHong Zhang bjtmp = bj + bi[i]; 17002efa7f71SHong Zhang batmp = ba + bs2*bi[i]; 17012efa7f71SHong Zhang /* apply level dropping rule to L part */ 17022efa7f71SHong Zhang ncut = nzi_al + dtcount; 17032efa7f71SHong Zhang if (ncut < nzi_bl) { 1704021822bcSHong Zhang ierr = PetscSortSplitReal(ncut,nzi_bl,vtmp_abs,jtmp);CHKERRQ(ierr); 17052efa7f71SHong Zhang ierr = PetscSortIntWithScalarArray(ncut,jtmp,vtmp);CHKERRQ(ierr); 17062efa7f71SHong Zhang } else { 17072efa7f71SHong Zhang ncut = nzi_bl; 17082efa7f71SHong Zhang } 17092efa7f71SHong Zhang for (j=0; j<ncut; j++) { 17102efa7f71SHong Zhang bjtmp[j] = jtmp[j]; 17112efa7f71SHong Zhang ierr = PetscMemcpy(batmp+bs2*j,rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 17122efa7f71SHong Zhang /* 17132efa7f71SHong Zhang printf(" col %d: ",bjtmp[j]); 17142efa7f71SHong Zhang for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*j+j1)); 17152efa7f71SHong Zhang printf("\n"); 17162efa7f71SHong Zhang */ 17172efa7f71SHong Zhang } 17182efa7f71SHong Zhang bi[i+1] = bi[i] + ncut; 17192efa7f71SHong Zhang nzi = ncut + 1; 17202efa7f71SHong Zhang 17212efa7f71SHong Zhang /* apply level dropping rule to U part */ 17222efa7f71SHong Zhang ncut = nzi_au + dtcount; 17232efa7f71SHong Zhang if (ncut < nzi_bu) { 1724021822bcSHong Zhang ierr = PetscSortSplitReal(ncut,nzi_bu,vtmp_abs+nzi_bl+1,jtmp+nzi_bl+1);CHKERRQ(ierr); 17252efa7f71SHong Zhang ierr = PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);CHKERRQ(ierr); 17262efa7f71SHong Zhang } else { 17272efa7f71SHong Zhang ncut = nzi_bu; 17282efa7f71SHong Zhang } 17292efa7f71SHong Zhang nzi += ncut; 17302efa7f71SHong Zhang 17312efa7f71SHong Zhang /* mark bdiagonal */ 17322efa7f71SHong Zhang bdiag[i+1] = bdiag[i] - (ncut + 1); 17336bce7ff8SHong Zhang bi[2*mbs - i] = bi[2*mbs - i +1] - (ncut + 1); 17346bce7ff8SHong Zhang 17352efa7f71SHong Zhang bjtmp = bj + bdiag[i]; 17362efa7f71SHong Zhang batmp = ba + bs2*bdiag[i]; 17372efa7f71SHong Zhang ierr = PetscMemcpy(batmp,rtmp+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 17382efa7f71SHong Zhang *bjtmp = i; 17392efa7f71SHong Zhang /* 17402efa7f71SHong Zhang printf(" diag %d: ",*bjtmp); 17412efa7f71SHong Zhang for (j=0; j<bs2; j++) { 17422efa7f71SHong Zhang printf(" %g,",batmp[j]); 17432efa7f71SHong Zhang } 17442efa7f71SHong Zhang printf("\n"); 17452efa7f71SHong Zhang */ 17462efa7f71SHong Zhang bjtmp = bj + bdiag[i+1]+1; 17472efa7f71SHong Zhang batmp = ba + (bdiag[i+1]+1)*bs2; 17482efa7f71SHong Zhang 17492efa7f71SHong Zhang for (k=0; k<ncut; k++) { 17502efa7f71SHong Zhang bjtmp[k] = jtmp[nzi_bl+1+k]; 17512efa7f71SHong Zhang ierr = PetscMemcpy(batmp+bs2*k,rtmp+bs2*bjtmp[k],bs2*sizeof(MatScalar));CHKERRQ(ierr); 17522efa7f71SHong Zhang /* 17532efa7f71SHong Zhang printf(" col %d:",bjtmp[k]); 17542efa7f71SHong Zhang for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*k+j1)); 17552efa7f71SHong Zhang printf("\n"); 17562efa7f71SHong Zhang */ 17572efa7f71SHong Zhang } 17582efa7f71SHong Zhang 17592efa7f71SHong Zhang im[i] = nzi; /* used by PetscLLAddSortedLU() */ 17602efa7f71SHong Zhang 17612efa7f71SHong Zhang /* invert diagonal block for simplier triangular solves - add shift??? */ 17622efa7f71SHong Zhang batmp = ba + bs2*bdiag[i]; 176396b95a6bSBarry Smith ierr = PetscKernel_A_gets_inverse_A(bs,batmp,v_pivots,v_work);CHKERRQ(ierr); 17642efa7f71SHong Zhang } /* for (i=0; i<mbs; i++) */ 1765fca92195SBarry Smith ierr = PetscFree3(v_work,multiplier,v_pivots);CHKERRQ(ierr); 17662efa7f71SHong Zhang 17676da515a1SHong Zhang /* printf("end of L %d, beginning of U %d\n",bi[mbs],bdiag[mbs]); */ 1768e32f2f54SBarry 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]); 17692efa7f71SHong Zhang 17702efa7f71SHong Zhang ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 17712efa7f71SHong Zhang ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 17722efa7f71SHong Zhang 17732efa7f71SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 17742efa7f71SHong Zhang 1775fca92195SBarry Smith ierr = PetscFree2(im,jtmp);CHKERRQ(ierr); 1776fca92195SBarry Smith ierr = PetscFree2(rtmp,vtmp);CHKERRQ(ierr); 17772efa7f71SHong Zhang 17782efa7f71SHong Zhang ierr = PetscLogFlops(bs2*B->cmap->n);CHKERRQ(ierr); 17792efa7f71SHong Zhang b->maxnz = b->nz = bi[mbs] + bdiag[0] - bdiag[mbs]; 17802efa7f71SHong Zhang 17812efa7f71SHong Zhang ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 17822efa7f71SHong Zhang ierr = ISIdentity(isicol,&icol_identity);CHKERRQ(ierr); 17832efa7f71SHong Zhang if (row_identity && icol_identity) { 17844dd39f65SShri Abhyankar B->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering; 17852efa7f71SHong Zhang } else { 17864dd39f65SShri Abhyankar B->ops->solve = MatSolve_SeqBAIJ_N; 17872efa7f71SHong Zhang } 17882efa7f71SHong Zhang 17892efa7f71SHong Zhang B->ops->solveadd = 0; 17902efa7f71SHong Zhang B->ops->solvetranspose = 0; 17912efa7f71SHong Zhang B->ops->solvetransposeadd = 0; 17922efa7f71SHong Zhang B->ops->matsolve = 0; 17932efa7f71SHong Zhang B->assembled = PETSC_TRUE; 17942efa7f71SHong Zhang B->preallocated = PETSC_TRUE; 1795c8342467SHong Zhang PetscFunctionReturn(0); 1796c8342467SHong Zhang } 1797