1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 383287d42SBarry Smith /* 483287d42SBarry Smith Factorization code for BAIJ format. 583287d42SBarry Smith */ 67c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h" 7c60f0209SBarry Smith #include "../src/mat/blockinvert.h" 883287d42SBarry Smith 983287d42SBarry Smith /* 1083287d42SBarry Smith Version for when blocks are 3 by 3 1183287d42SBarry Smith */ 124a2ae208SSatish Balay #undef __FUNCT__ 1306e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_inplace" 1406e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_inplace(Mat C,Mat A,const MatFactorInfo *info) 1583287d42SBarry Smith { 1683287d42SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 1783287d42SBarry Smith IS isrow = b->row,isicol = b->icol; 186849ba73SBarry Smith PetscErrorCode ierr; 195d0c19d7SBarry Smith const PetscInt *r,*ic; 205d0c19d7SBarry Smith PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 21690b6cddSBarry Smith PetscInt *ajtmpold,*ajtmp,nz,row,*ai=a->i,*aj=a->j; 22690b6cddSBarry Smith PetscInt *diag_offset = b->diag,idx,*pj; 2383287d42SBarry Smith MatScalar *pv,*v,*rtmp,*pc,*w,*x; 2483287d42SBarry Smith MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4; 2583287d42SBarry Smith MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9; 2683287d42SBarry Smith MatScalar *ba = b->a,*aa = a->a; 2762bba022SBarry Smith PetscReal shift = info->shiftinblocks; 2883287d42SBarry Smith 2983287d42SBarry Smith PetscFunctionBegin; 3083287d42SBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 3183287d42SBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 32b0a32e0cSBarry Smith ierr = PetscMalloc(9*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 3383287d42SBarry Smith 3483287d42SBarry Smith for (i=0; i<n; i++) { 3583287d42SBarry Smith nz = bi[i+1] - bi[i]; 3683287d42SBarry Smith ajtmp = bj + bi[i]; 3783287d42SBarry Smith for (j=0; j<nz; j++) { 3883287d42SBarry Smith x = rtmp + 9*ajtmp[j]; 3983287d42SBarry Smith x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = 0.0; 4083287d42SBarry Smith } 4183287d42SBarry Smith /* load in initial (unfactored row) */ 4283287d42SBarry Smith idx = r[i]; 4383287d42SBarry Smith nz = ai[idx+1] - ai[idx]; 4483287d42SBarry Smith ajtmpold = aj + ai[idx]; 4583287d42SBarry Smith v = aa + 9*ai[idx]; 4683287d42SBarry Smith for (j=0; j<nz; j++) { 4783287d42SBarry Smith x = rtmp + 9*ic[ajtmpold[j]]; 4883287d42SBarry Smith x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 4983287d42SBarry Smith x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8]; 5083287d42SBarry Smith v += 9; 5183287d42SBarry Smith } 5283287d42SBarry Smith row = *ajtmp++; 5383287d42SBarry Smith while (row < i) { 5483287d42SBarry Smith pc = rtmp + 9*row; 5583287d42SBarry Smith p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 5683287d42SBarry Smith p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8]; 5783287d42SBarry Smith if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || 5883287d42SBarry Smith p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) { 5983287d42SBarry Smith pv = ba + 9*diag_offset[row]; 6083287d42SBarry Smith pj = bj + diag_offset[row] + 1; 6183287d42SBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 6283287d42SBarry Smith x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 6383287d42SBarry Smith pc[0] = m1 = p1*x1 + p4*x2 + p7*x3; 6483287d42SBarry Smith pc[1] = m2 = p2*x1 + p5*x2 + p8*x3; 6583287d42SBarry Smith pc[2] = m3 = p3*x1 + p6*x2 + p9*x3; 6683287d42SBarry Smith 6783287d42SBarry Smith pc[3] = m4 = p1*x4 + p4*x5 + p7*x6; 6883287d42SBarry Smith pc[4] = m5 = p2*x4 + p5*x5 + p8*x6; 6983287d42SBarry Smith pc[5] = m6 = p3*x4 + p6*x5 + p9*x6; 7083287d42SBarry Smith 7183287d42SBarry Smith pc[6] = m7 = p1*x7 + p4*x8 + p7*x9; 7283287d42SBarry Smith pc[7] = m8 = p2*x7 + p5*x8 + p8*x9; 7383287d42SBarry Smith pc[8] = m9 = p3*x7 + p6*x8 + p9*x9; 7483287d42SBarry Smith nz = bi[row+1] - diag_offset[row] - 1; 7583287d42SBarry Smith pv += 9; 7683287d42SBarry Smith for (j=0; j<nz; j++) { 7783287d42SBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 7883287d42SBarry Smith x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 7983287d42SBarry Smith x = rtmp + 9*pj[j]; 8083287d42SBarry Smith x[0] -= m1*x1 + m4*x2 + m7*x3; 8183287d42SBarry Smith x[1] -= m2*x1 + m5*x2 + m8*x3; 8283287d42SBarry Smith x[2] -= m3*x1 + m6*x2 + m9*x3; 8383287d42SBarry Smith 8483287d42SBarry Smith x[3] -= m1*x4 + m4*x5 + m7*x6; 8583287d42SBarry Smith x[4] -= m2*x4 + m5*x5 + m8*x6; 8683287d42SBarry Smith x[5] -= m3*x4 + m6*x5 + m9*x6; 8783287d42SBarry Smith 8883287d42SBarry Smith x[6] -= m1*x7 + m4*x8 + m7*x9; 8983287d42SBarry Smith x[7] -= m2*x7 + m5*x8 + m8*x9; 9083287d42SBarry Smith x[8] -= m3*x7 + m6*x8 + m9*x9; 9183287d42SBarry Smith pv += 9; 9283287d42SBarry Smith } 93dc0b31edSSatish Balay ierr = PetscLogFlops(54.0*nz+36.0);CHKERRQ(ierr); 9483287d42SBarry Smith } 9583287d42SBarry Smith row = *ajtmp++; 9683287d42SBarry Smith } 9783287d42SBarry Smith /* finished row so stick it into b->a */ 9883287d42SBarry Smith pv = ba + 9*bi[i]; 9983287d42SBarry Smith pj = bj + bi[i]; 10083287d42SBarry Smith nz = bi[i+1] - bi[i]; 10183287d42SBarry Smith for (j=0; j<nz; j++) { 10283287d42SBarry Smith x = rtmp + 9*pj[j]; 10383287d42SBarry Smith pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 10483287d42SBarry Smith pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8]; 10583287d42SBarry Smith pv += 9; 10683287d42SBarry Smith } 10783287d42SBarry Smith /* invert diagonal block */ 10883287d42SBarry Smith w = ba + 9*diag_offset[i]; 10962bba022SBarry Smith ierr = Kernel_A_gets_inverse_A_3(w,shift);CHKERRQ(ierr); 11083287d42SBarry Smith } 11183287d42SBarry Smith 11283287d42SBarry Smith ierr = PetscFree(rtmp);CHKERRQ(ierr); 11383287d42SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 11483287d42SBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 11506e38f1dSHong Zhang C->ops->solve = MatSolve_SeqBAIJ_3_inplace; 11606e38f1dSHong Zhang C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_inplace; 11783287d42SBarry Smith C->assembled = PETSC_TRUE; 118766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*3*3*3*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 11983287d42SBarry Smith PetscFunctionReturn(0); 12083287d42SBarry Smith } 12117542729SShri Abhyankar 1224dd39f65SShri Abhyankar /* MatLUFactorNumeric_SeqBAIJ_3 - 1234dd39f65SShri Abhyankar copied from MatLUFactorNumeric_SeqBAIJ_N_inplace() and manually re-implemented 12417542729SShri Abhyankar Kernel_A_gets_A_times_B() 12517542729SShri Abhyankar Kernel_A_gets_A_minus_B_times_C() 12617542729SShri Abhyankar Kernel_A_gets_inverse_A() 12717542729SShri Abhyankar */ 12817542729SShri Abhyankar #undef __FUNCT__ 1294dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3" 1304dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3(Mat B,Mat A,const MatFactorInfo *info) 13117542729SShri Abhyankar { 13217542729SShri Abhyankar Mat C=B; 13317542729SShri Abhyankar Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; 13417542729SShri Abhyankar IS isrow = b->row,isicol = b->icol; 13517542729SShri Abhyankar PetscErrorCode ierr; 13617542729SShri Abhyankar const PetscInt *r,*ic,*ics; 137*bbd65245SShri Abhyankar PetscInt i,j,k,nz,nzL,row; 138*bbd65245SShri Abhyankar const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 139*bbd65245SShri Abhyankar const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2; 14017542729SShri Abhyankar MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 141*bbd65245SShri Abhyankar PetscInt flg; 14217542729SShri Abhyankar PetscReal shift = info->shiftinblocks; 14317542729SShri Abhyankar 14417542729SShri Abhyankar PetscFunctionBegin; 14517542729SShri Abhyankar ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 14617542729SShri Abhyankar ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 14717542729SShri Abhyankar 14817542729SShri Abhyankar /* generate work space needed by the factorization */ 149fca92195SBarry Smith ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr); 15017542729SShri Abhyankar ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 15117542729SShri Abhyankar ics = ic; 15217542729SShri Abhyankar 15317542729SShri Abhyankar for (i=0; i<n; i++){ 15417542729SShri Abhyankar /* zero rtmp */ 15517542729SShri Abhyankar /* L part */ 15617542729SShri Abhyankar nz = bi[i+1] - bi[i]; 15717542729SShri Abhyankar bjtmp = bj + bi[i]; 15817542729SShri Abhyankar for (j=0; j<nz; j++){ 15917542729SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 16017542729SShri Abhyankar } 16117542729SShri Abhyankar 16217542729SShri Abhyankar /* U part */ 1630c4413a7SShri Abhyankar nz = bdiag[i] - bdiag[i+1]; 1640c4413a7SShri Abhyankar bjtmp = bj + bdiag[i+1]+1; 1650c4413a7SShri Abhyankar for (j=0; j<nz; j++){ 1660c4413a7SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1670c4413a7SShri Abhyankar } 1680c4413a7SShri Abhyankar 1690c4413a7SShri Abhyankar /* load in initial (unfactored row) */ 1700c4413a7SShri Abhyankar nz = ai[r[i]+1] - ai[r[i]]; 1710c4413a7SShri Abhyankar ajtmp = aj + ai[r[i]]; 1720c4413a7SShri Abhyankar v = aa + bs2*ai[r[i]]; 1730c4413a7SShri Abhyankar for (j=0; j<nz; j++) { 1740c4413a7SShri Abhyankar ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1750c4413a7SShri Abhyankar } 1760c4413a7SShri Abhyankar 1770c4413a7SShri Abhyankar /* elimination */ 1780c4413a7SShri Abhyankar bjtmp = bj + bi[i]; 1790c4413a7SShri Abhyankar nzL = bi[i+1] - bi[i]; 1800c4413a7SShri Abhyankar for(k = 0;k < nzL;k++){ 1810c4413a7SShri Abhyankar row = bjtmp[k]; 1820c4413a7SShri Abhyankar pc = rtmp + bs2*row; 1830c4413a7SShri Abhyankar for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }} 1840c4413a7SShri Abhyankar if (flg) { 1850c4413a7SShri Abhyankar pv = b->a + bs2*bdiag[row]; 1860c4413a7SShri Abhyankar /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 1870c4413a7SShri Abhyankar ierr = Kernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr); 1880c4413a7SShri Abhyankar 1890c4413a7SShri Abhyankar pj = b->j + bdiag[row+1] + 1; /* beginning of U(row,:) */ 1900c4413a7SShri Abhyankar pv = b->a + bs2*(bdiag[row+1]+1); 1910c4413a7SShri Abhyankar nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */ 1920c4413a7SShri Abhyankar for (j=0; j<nz; j++) { 1930c4413a7SShri Abhyankar /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 1940c4413a7SShri Abhyankar /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 1950c4413a7SShri Abhyankar v = rtmp + bs2*pj[j]; 1960c4413a7SShri Abhyankar ierr = Kernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr); 1970c4413a7SShri Abhyankar pv += bs2; 1980c4413a7SShri Abhyankar } 1990c4413a7SShri Abhyankar ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 2000c4413a7SShri Abhyankar } 2010c4413a7SShri Abhyankar } 2020c4413a7SShri Abhyankar 2030c4413a7SShri Abhyankar /* finished row so stick it into b->a */ 2040c4413a7SShri Abhyankar /* L part */ 2050c4413a7SShri Abhyankar pv = b->a + bs2*bi[i] ; 2060c4413a7SShri Abhyankar pj = b->j + bi[i] ; 2070c4413a7SShri Abhyankar nz = bi[i+1] - bi[i]; 2080c4413a7SShri Abhyankar for (j=0; j<nz; j++) { 2090c4413a7SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 2100c4413a7SShri Abhyankar } 2110c4413a7SShri Abhyankar 2120c4413a7SShri Abhyankar /* Mark diagonal and invert diagonal for simplier triangular solves */ 2130c4413a7SShri Abhyankar pv = b->a + bs2*bdiag[i]; 2140c4413a7SShri Abhyankar pj = b->j + bdiag[i]; 2150c4413a7SShri Abhyankar ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 2160c4413a7SShri Abhyankar /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 2170c4413a7SShri Abhyankar ierr = Kernel_A_gets_inverse_A_3(pv,shift);CHKERRQ(ierr); 2180c4413a7SShri Abhyankar 2190c4413a7SShri Abhyankar /* U part */ 2200c4413a7SShri Abhyankar pj = b->j + bdiag[i+1] + 1; 2210c4413a7SShri Abhyankar pv = b->a + bs2*(bdiag[i+1]+1); 2220c4413a7SShri Abhyankar nz = bdiag[i] - bdiag[i+1] - 1; 2230c4413a7SShri Abhyankar for (j=0; j<nz; j++){ 2240c4413a7SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 2250c4413a7SShri Abhyankar } 2260c4413a7SShri Abhyankar } 2270c4413a7SShri Abhyankar 228fca92195SBarry Smith ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); 2290c4413a7SShri Abhyankar ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 2300c4413a7SShri Abhyankar ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 2314dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_3; 2324dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3; 2330c4413a7SShri Abhyankar 2340c4413a7SShri Abhyankar C->assembled = PETSC_TRUE; 235766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*3*3*3*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 2360c4413a7SShri Abhyankar PetscFunctionReturn(0); 2370c4413a7SShri Abhyankar } 2380c4413a7SShri Abhyankar 2390c4413a7SShri Abhyankar #undef __FUNCT__ 24006e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace" 24106e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info) 24217542729SShri Abhyankar { 24317542729SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 24417542729SShri Abhyankar PetscErrorCode ierr; 24517542729SShri Abhyankar PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 24617542729SShri Abhyankar PetscInt *ajtmpold,*ajtmp,nz,row; 24717542729SShri Abhyankar PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; 24817542729SShri Abhyankar MatScalar *pv,*v,*rtmp,*pc,*w,*x; 24917542729SShri Abhyankar MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4; 25017542729SShri Abhyankar MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9; 25117542729SShri Abhyankar MatScalar *ba = b->a,*aa = a->a; 25217542729SShri Abhyankar PetscReal shift = info->shiftinblocks; 25317542729SShri Abhyankar 25417542729SShri Abhyankar PetscFunctionBegin; 25517542729SShri Abhyankar ierr = PetscMalloc(9*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 25617542729SShri Abhyankar 25717542729SShri Abhyankar for (i=0; i<n; i++) { 25817542729SShri Abhyankar nz = bi[i+1] - bi[i]; 25917542729SShri Abhyankar ajtmp = bj + bi[i]; 26017542729SShri Abhyankar for (j=0; j<nz; j++) { 26117542729SShri Abhyankar x = rtmp+9*ajtmp[j]; 26217542729SShri Abhyankar x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = 0.0; 26317542729SShri Abhyankar } 26417542729SShri Abhyankar /* load in initial (unfactored row) */ 26517542729SShri Abhyankar nz = ai[i+1] - ai[i]; 26617542729SShri Abhyankar ajtmpold = aj + ai[i]; 26717542729SShri Abhyankar v = aa + 9*ai[i]; 26817542729SShri Abhyankar for (j=0; j<nz; j++) { 26917542729SShri Abhyankar x = rtmp+9*ajtmpold[j]; 27017542729SShri Abhyankar x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 27117542729SShri Abhyankar x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8]; 27217542729SShri Abhyankar v += 9; 27317542729SShri Abhyankar } 27417542729SShri Abhyankar row = *ajtmp++; 27517542729SShri Abhyankar while (row < i) { 27617542729SShri Abhyankar pc = rtmp + 9*row; 27717542729SShri Abhyankar p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 27817542729SShri Abhyankar p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8]; 27917542729SShri Abhyankar if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || 28017542729SShri Abhyankar p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) { 28117542729SShri Abhyankar pv = ba + 9*diag_offset[row]; 28217542729SShri Abhyankar pj = bj + diag_offset[row] + 1; 28317542729SShri Abhyankar x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 28417542729SShri Abhyankar x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 28517542729SShri Abhyankar pc[0] = m1 = p1*x1 + p4*x2 + p7*x3; 28617542729SShri Abhyankar pc[1] = m2 = p2*x1 + p5*x2 + p8*x3; 28717542729SShri Abhyankar pc[2] = m3 = p3*x1 + p6*x2 + p9*x3; 28817542729SShri Abhyankar 28917542729SShri Abhyankar pc[3] = m4 = p1*x4 + p4*x5 + p7*x6; 29017542729SShri Abhyankar pc[4] = m5 = p2*x4 + p5*x5 + p8*x6; 29117542729SShri Abhyankar pc[5] = m6 = p3*x4 + p6*x5 + p9*x6; 29217542729SShri Abhyankar 29317542729SShri Abhyankar pc[6] = m7 = p1*x7 + p4*x8 + p7*x9; 29417542729SShri Abhyankar pc[7] = m8 = p2*x7 + p5*x8 + p8*x9; 29517542729SShri Abhyankar pc[8] = m9 = p3*x7 + p6*x8 + p9*x9; 29617542729SShri Abhyankar 29717542729SShri Abhyankar nz = bi[row+1] - diag_offset[row] - 1; 29817542729SShri Abhyankar pv += 9; 29917542729SShri Abhyankar for (j=0; j<nz; j++) { 30017542729SShri Abhyankar x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 30117542729SShri Abhyankar x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 30217542729SShri Abhyankar x = rtmp + 9*pj[j]; 30317542729SShri Abhyankar x[0] -= m1*x1 + m4*x2 + m7*x3; 30417542729SShri Abhyankar x[1] -= m2*x1 + m5*x2 + m8*x3; 30517542729SShri Abhyankar x[2] -= m3*x1 + m6*x2 + m9*x3; 30617542729SShri Abhyankar 30717542729SShri Abhyankar x[3] -= m1*x4 + m4*x5 + m7*x6; 30817542729SShri Abhyankar x[4] -= m2*x4 + m5*x5 + m8*x6; 30917542729SShri Abhyankar x[5] -= m3*x4 + m6*x5 + m9*x6; 31017542729SShri Abhyankar 31117542729SShri Abhyankar x[6] -= m1*x7 + m4*x8 + m7*x9; 31217542729SShri Abhyankar x[7] -= m2*x7 + m5*x8 + m8*x9; 31317542729SShri Abhyankar x[8] -= m3*x7 + m6*x8 + m9*x9; 31417542729SShri Abhyankar pv += 9; 31517542729SShri Abhyankar } 31617542729SShri Abhyankar ierr = PetscLogFlops(54.0*nz+36.0);CHKERRQ(ierr); 31717542729SShri Abhyankar } 31817542729SShri Abhyankar row = *ajtmp++; 31917542729SShri Abhyankar } 32017542729SShri Abhyankar /* finished row so stick it into b->a */ 32117542729SShri Abhyankar pv = ba + 9*bi[i]; 32217542729SShri Abhyankar pj = bj + bi[i]; 32317542729SShri Abhyankar nz = bi[i+1] - bi[i]; 32417542729SShri Abhyankar for (j=0; j<nz; j++) { 32517542729SShri Abhyankar x = rtmp+9*pj[j]; 32617542729SShri Abhyankar pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 32717542729SShri Abhyankar pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8]; 32817542729SShri Abhyankar pv += 9; 32917542729SShri Abhyankar } 33017542729SShri Abhyankar /* invert diagonal block */ 33117542729SShri Abhyankar w = ba + 9*diag_offset[i]; 33217542729SShri Abhyankar ierr = Kernel_A_gets_inverse_A_3(w,shift);CHKERRQ(ierr); 33317542729SShri Abhyankar } 33417542729SShri Abhyankar 33517542729SShri Abhyankar ierr = PetscFree(rtmp);CHKERRQ(ierr); 33606e38f1dSHong Zhang C->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering_inplace; 33706e38f1dSHong Zhang C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering_inplace; 33817542729SShri Abhyankar C->assembled = PETSC_TRUE; 339766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*3*3*3*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 34017542729SShri Abhyankar PetscFunctionReturn(0); 34117542729SShri Abhyankar } 34217542729SShri Abhyankar 34317542729SShri Abhyankar /* 3444dd39f65SShri Abhyankar MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering - 3454dd39f65SShri Abhyankar copied from MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace() 34617542729SShri Abhyankar */ 34717542729SShri Abhyankar #undef __FUNCT__ 3484dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering" 3494dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info) 35017542729SShri Abhyankar { 35117542729SShri Abhyankar Mat C=B; 35217542729SShri Abhyankar Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; 35317542729SShri Abhyankar PetscErrorCode ierr; 354*bbd65245SShri Abhyankar PetscInt i,j,k,nz,nzL,row; 355*bbd65245SShri Abhyankar const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 356*bbd65245SShri Abhyankar const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2; 35717542729SShri Abhyankar MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 358*bbd65245SShri Abhyankar PetscInt flg; 35917542729SShri Abhyankar PetscReal shift = info->shiftinblocks; 36017542729SShri Abhyankar 36117542729SShri Abhyankar PetscFunctionBegin; 36217542729SShri Abhyankar /* generate work space needed by the factorization */ 363fca92195SBarry Smith ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr); 36417542729SShri Abhyankar ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 36517542729SShri Abhyankar 36617542729SShri Abhyankar for (i=0; i<n; i++){ 36717542729SShri Abhyankar /* zero rtmp */ 36817542729SShri Abhyankar /* L part */ 36917542729SShri Abhyankar nz = bi[i+1] - bi[i]; 37017542729SShri Abhyankar bjtmp = bj + bi[i]; 37117542729SShri Abhyankar for (j=0; j<nz; j++){ 37217542729SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 37317542729SShri Abhyankar } 37417542729SShri Abhyankar 37517542729SShri Abhyankar /* U part */ 376b2b2dd24SShri Abhyankar nz = bdiag[i] - bdiag[i+1]; 377b2b2dd24SShri Abhyankar bjtmp = bj + bdiag[i+1] + 1; 378b2b2dd24SShri Abhyankar for (j=0; j<nz; j++){ 379b2b2dd24SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 380b2b2dd24SShri Abhyankar } 381b2b2dd24SShri Abhyankar 382b2b2dd24SShri Abhyankar /* load in initial (unfactored row) */ 383b2b2dd24SShri Abhyankar nz = ai[i+1] - ai[i]; 384b2b2dd24SShri Abhyankar ajtmp = aj + ai[i]; 385b2b2dd24SShri Abhyankar v = aa + bs2*ai[i]; 386b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 387b2b2dd24SShri Abhyankar ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 388b2b2dd24SShri Abhyankar } 389b2b2dd24SShri Abhyankar 390b2b2dd24SShri Abhyankar /* elimination */ 391b2b2dd24SShri Abhyankar bjtmp = bj + bi[i]; 392b2b2dd24SShri Abhyankar nzL = bi[i+1] - bi[i]; 393b2b2dd24SShri Abhyankar for(k=0;k<nzL;k++){ 394b2b2dd24SShri Abhyankar row = bjtmp[k]; 395b2b2dd24SShri Abhyankar pc = rtmp + bs2*row; 396b2b2dd24SShri Abhyankar for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }} 397b2b2dd24SShri Abhyankar if (flg) { 398b2b2dd24SShri Abhyankar pv = b->a + bs2*bdiag[row]; 399b2b2dd24SShri Abhyankar /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 400b2b2dd24SShri Abhyankar ierr = Kernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr); 401b2b2dd24SShri Abhyankar 402b2b2dd24SShri Abhyankar pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */ 403b2b2dd24SShri Abhyankar pv = b->a + bs2*(bdiag[row+1]+1); 404b2b2dd24SShri Abhyankar nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */ 405b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 406b2b2dd24SShri Abhyankar /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 407b2b2dd24SShri Abhyankar /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 408b2b2dd24SShri Abhyankar v = rtmp + bs2*pj[j]; 409b2b2dd24SShri Abhyankar ierr = Kernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr); 410b2b2dd24SShri Abhyankar pv += bs2; 411b2b2dd24SShri Abhyankar } 412b2b2dd24SShri Abhyankar ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 413b2b2dd24SShri Abhyankar } 414b2b2dd24SShri Abhyankar } 415b2b2dd24SShri Abhyankar 416b2b2dd24SShri Abhyankar /* finished row so stick it into b->a */ 417b2b2dd24SShri Abhyankar /* L part */ 418b2b2dd24SShri Abhyankar pv = b->a + bs2*bi[i] ; 419b2b2dd24SShri Abhyankar pj = b->j + bi[i] ; 420b2b2dd24SShri Abhyankar nz = bi[i+1] - bi[i]; 421b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 422b2b2dd24SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 423b2b2dd24SShri Abhyankar } 424b2b2dd24SShri Abhyankar 425b2b2dd24SShri Abhyankar /* Mark diagonal and invert diagonal for simplier triangular solves */ 426b2b2dd24SShri Abhyankar pv = b->a + bs2*bdiag[i]; 427b2b2dd24SShri Abhyankar pj = b->j + bdiag[i]; 428b2b2dd24SShri Abhyankar ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 429b2b2dd24SShri Abhyankar /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 430b2b2dd24SShri Abhyankar ierr = Kernel_A_gets_inverse_A_3(pv,shift);CHKERRQ(ierr); 431b2b2dd24SShri Abhyankar 432b2b2dd24SShri Abhyankar /* U part */ 433b2b2dd24SShri Abhyankar pv = b->a + bs2*(bdiag[i+1]+1); 434b2b2dd24SShri Abhyankar pj = b->j + bdiag[i+1]+1; 435b2b2dd24SShri Abhyankar nz = bdiag[i] - bdiag[i+1] - 1; 436b2b2dd24SShri Abhyankar for (j=0; j<nz; j++){ 437b2b2dd24SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 438b2b2dd24SShri Abhyankar } 439b2b2dd24SShri Abhyankar } 440fca92195SBarry Smith ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); 4414dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering; 4424dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering; 443b2b2dd24SShri Abhyankar C->assembled = PETSC_TRUE; 444766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*3*3*3*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 445b2b2dd24SShri Abhyankar PetscFunctionReturn(0); 446b2b2dd24SShri Abhyankar } 447b2b2dd24SShri Abhyankar 448