1be1d678aSKris Buschelman 283287d42SBarry Smith /* 383287d42SBarry Smith Factorization code for BAIJ format. 483287d42SBarry Smith */ 5c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h> 6af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h> 783287d42SBarry Smith 883287d42SBarry Smith /* 983287d42SBarry Smith Version for when blocks are 3 by 3 1083287d42SBarry Smith */ 114a2ae208SSatish Balay #undef __FUNCT__ 1206e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_inplace" 1306e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_inplace(Mat C,Mat A,const MatFactorInfo *info) 1483287d42SBarry Smith { 1583287d42SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data; 1683287d42SBarry Smith IS isrow = b->row,isicol = b->icol; 176849ba73SBarry Smith PetscErrorCode ierr; 185d0c19d7SBarry Smith const PetscInt *r,*ic; 195d0c19d7SBarry Smith PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 20690b6cddSBarry Smith PetscInt *ajtmpold,*ajtmp,nz,row,*ai=a->i,*aj=a->j; 21690b6cddSBarry Smith PetscInt *diag_offset = b->diag,idx,*pj; 2283287d42SBarry Smith MatScalar *pv,*v,*rtmp,*pc,*w,*x; 2383287d42SBarry Smith MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4; 2483287d42SBarry Smith MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9; 2583287d42SBarry Smith MatScalar *ba = b->a,*aa = a->a; 26182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 27*a455e926SHong Zhang PetscBool allowzeropivot,zeropivotdetected; 2883287d42SBarry Smith 2983287d42SBarry Smith PetscFunctionBegin; 3083287d42SBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 3183287d42SBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 32785e854fSJed Brown ierr = PetscMalloc1(9*(n+1),&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]; 109*a455e926SHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 110*a455e926SHong Zhang ierr = PetscKernel_A_gets_inverse_A_3(w,shift,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr); 111922032d7SHong Zhang if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 11283287d42SBarry Smith } 11383287d42SBarry Smith 11483287d42SBarry Smith ierr = PetscFree(rtmp);CHKERRQ(ierr); 11583287d42SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 11683287d42SBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 11726fbe8dcSKarl Rupp 11806e38f1dSHong Zhang C->ops->solve = MatSolve_SeqBAIJ_3_inplace; 11906e38f1dSHong Zhang C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_inplace; 12083287d42SBarry Smith C->assembled = PETSC_TRUE; 12126fbe8dcSKarl Rupp 122766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*3*3*3*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 12383287d42SBarry Smith PetscFunctionReturn(0); 12483287d42SBarry Smith } 12517542729SShri Abhyankar 1264dd39f65SShri Abhyankar /* MatLUFactorNumeric_SeqBAIJ_3 - 1274dd39f65SShri Abhyankar copied from MatLUFactorNumeric_SeqBAIJ_N_inplace() and manually re-implemented 12896b95a6bSBarry Smith PetscKernel_A_gets_A_times_B() 12996b95a6bSBarry Smith PetscKernel_A_gets_A_minus_B_times_C() 13096b95a6bSBarry Smith PetscKernel_A_gets_inverse_A() 13117542729SShri Abhyankar */ 13217542729SShri Abhyankar #undef __FUNCT__ 1334dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3" 1344dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3(Mat B,Mat A,const MatFactorInfo *info) 13517542729SShri Abhyankar { 13617542729SShri Abhyankar Mat C =B; 13717542729SShri Abhyankar Mat_SeqBAIJ *a =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data; 13817542729SShri Abhyankar IS isrow = b->row,isicol = b->icol; 13917542729SShri Abhyankar PetscErrorCode ierr; 1405a586d82SBarry Smith const PetscInt *r,*ic; 141bbd65245SShri Abhyankar PetscInt i,j,k,nz,nzL,row; 142bbd65245SShri Abhyankar const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 143bbd65245SShri Abhyankar const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2; 14417542729SShri Abhyankar MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 145bbd65245SShri Abhyankar PetscInt flg; 146182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 147*a455e926SHong Zhang PetscBool allowzeropivot,zeropivotdetected; 14817542729SShri Abhyankar 14917542729SShri Abhyankar PetscFunctionBegin; 15017542729SShri Abhyankar ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 15117542729SShri Abhyankar ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 15217542729SShri Abhyankar 15317542729SShri Abhyankar /* generate work space needed by the factorization */ 154dcca6d9dSJed Brown ierr = PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);CHKERRQ(ierr); 15517542729SShri Abhyankar ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 15617542729SShri Abhyankar 15717542729SShri Abhyankar for (i=0; i<n; i++) { 15817542729SShri Abhyankar /* zero rtmp */ 15917542729SShri Abhyankar /* L part */ 16017542729SShri Abhyankar nz = bi[i+1] - bi[i]; 16117542729SShri Abhyankar bjtmp = bj + bi[i]; 16217542729SShri Abhyankar for (j=0; j<nz; j++) { 16317542729SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 16417542729SShri Abhyankar } 16517542729SShri Abhyankar 16617542729SShri Abhyankar /* U part */ 1670c4413a7SShri Abhyankar nz = bdiag[i] - bdiag[i+1]; 1680c4413a7SShri Abhyankar bjtmp = bj + bdiag[i+1]+1; 1690c4413a7SShri Abhyankar for (j=0; j<nz; j++) { 1700c4413a7SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1710c4413a7SShri Abhyankar } 1720c4413a7SShri Abhyankar 1730c4413a7SShri Abhyankar /* load in initial (unfactored row) */ 1740c4413a7SShri Abhyankar nz = ai[r[i]+1] - ai[r[i]]; 1750c4413a7SShri Abhyankar ajtmp = aj + ai[r[i]]; 1760c4413a7SShri Abhyankar v = aa + bs2*ai[r[i]]; 1770c4413a7SShri Abhyankar for (j=0; j<nz; j++) { 1780c4413a7SShri Abhyankar ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1790c4413a7SShri Abhyankar } 1800c4413a7SShri Abhyankar 1810c4413a7SShri Abhyankar /* elimination */ 1820c4413a7SShri Abhyankar bjtmp = bj + bi[i]; 1830c4413a7SShri Abhyankar nzL = bi[i+1] - bi[i]; 1840c4413a7SShri Abhyankar for (k = 0; k < nzL; k++) { 1850c4413a7SShri Abhyankar row = bjtmp[k]; 1860c4413a7SShri Abhyankar pc = rtmp + bs2*row; 187c35f09e5SBarry Smith for (flg=0,j=0; j<bs2; j++) { 188c35f09e5SBarry Smith if (pc[j]!=0.0) { 189c35f09e5SBarry Smith flg = 1; 190c35f09e5SBarry Smith break; 191c35f09e5SBarry Smith } 192c35f09e5SBarry Smith } 1930c4413a7SShri Abhyankar if (flg) { 1940c4413a7SShri Abhyankar pv = b->a + bs2*bdiag[row]; 19596b95a6bSBarry Smith /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 19696b95a6bSBarry Smith ierr = PetscKernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr); 1970c4413a7SShri Abhyankar 1980c4413a7SShri Abhyankar pj = b->j + bdiag[row+1] + 1; /* beginning of U(row,:) */ 1990c4413a7SShri Abhyankar pv = b->a + bs2*(bdiag[row+1]+1); 2000c4413a7SShri Abhyankar nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */ 2010c4413a7SShri Abhyankar for (j=0; j<nz; j++) { 20296b95a6bSBarry Smith /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 2030c4413a7SShri Abhyankar /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 2040c4413a7SShri Abhyankar v = rtmp + bs2*pj[j]; 20596b95a6bSBarry Smith ierr = PetscKernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr); 2060c4413a7SShri Abhyankar pv += bs2; 2070c4413a7SShri Abhyankar } 2080c4413a7SShri Abhyankar ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 2090c4413a7SShri Abhyankar } 2100c4413a7SShri Abhyankar } 2110c4413a7SShri Abhyankar 2120c4413a7SShri Abhyankar /* finished row so stick it into b->a */ 2130c4413a7SShri Abhyankar /* L part */ 2140c4413a7SShri Abhyankar pv = b->a + bs2*bi[i]; 2150c4413a7SShri Abhyankar pj = b->j + bi[i]; 2160c4413a7SShri Abhyankar nz = bi[i+1] - bi[i]; 2170c4413a7SShri Abhyankar for (j=0; j<nz; j++) { 2180c4413a7SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 2190c4413a7SShri Abhyankar } 2200c4413a7SShri Abhyankar 2210c4413a7SShri Abhyankar /* Mark diagonal and invert diagonal for simplier triangular solves */ 2220c4413a7SShri Abhyankar pv = b->a + bs2*bdiag[i]; 2230c4413a7SShri Abhyankar pj = b->j + bdiag[i]; 2240c4413a7SShri Abhyankar ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 22596b95a6bSBarry Smith /* ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 226*a455e926SHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 227*a455e926SHong Zhang ierr = PetscKernel_A_gets_inverse_A_3(pv,shift,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr); 2282e92ee13SHong Zhang if (zeropivotdetected) B->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 2290c4413a7SShri Abhyankar 2300c4413a7SShri Abhyankar /* U part */ 2310c4413a7SShri Abhyankar pj = b->j + bdiag[i+1] + 1; 2320c4413a7SShri Abhyankar pv = b->a + bs2*(bdiag[i+1]+1); 2330c4413a7SShri Abhyankar nz = bdiag[i] - bdiag[i+1] - 1; 2340c4413a7SShri Abhyankar for (j=0; j<nz; j++) { 2350c4413a7SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 2360c4413a7SShri Abhyankar } 2370c4413a7SShri Abhyankar } 2380c4413a7SShri Abhyankar 239fca92195SBarry Smith ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); 2400c4413a7SShri Abhyankar ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 2410c4413a7SShri Abhyankar ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 24226fbe8dcSKarl Rupp 2434dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_3; 2444dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3; 2450c4413a7SShri Abhyankar C->assembled = PETSC_TRUE; 24626fbe8dcSKarl Rupp 247766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*3*3*3*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 2480c4413a7SShri Abhyankar PetscFunctionReturn(0); 2490c4413a7SShri Abhyankar } 2500c4413a7SShri Abhyankar 2510c4413a7SShri Abhyankar #undef __FUNCT__ 25206e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace" 25306e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info) 25417542729SShri Abhyankar { 25517542729SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data; 25617542729SShri Abhyankar PetscErrorCode ierr; 25717542729SShri Abhyankar PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 25817542729SShri Abhyankar PetscInt *ajtmpold,*ajtmp,nz,row; 25917542729SShri Abhyankar PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; 26017542729SShri Abhyankar MatScalar *pv,*v,*rtmp,*pc,*w,*x; 26117542729SShri Abhyankar MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4; 26217542729SShri Abhyankar MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9; 26317542729SShri Abhyankar MatScalar *ba = b->a,*aa = a->a; 264182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 265*a455e926SHong Zhang PetscBool allowzeropivot,zeropivotdetected; 26617542729SShri Abhyankar 26717542729SShri Abhyankar PetscFunctionBegin; 268785e854fSJed Brown ierr = PetscMalloc1(9*(n+1),&rtmp);CHKERRQ(ierr); 26917542729SShri Abhyankar 27017542729SShri Abhyankar for (i=0; i<n; i++) { 27117542729SShri Abhyankar nz = bi[i+1] - bi[i]; 27217542729SShri Abhyankar ajtmp = bj + bi[i]; 27317542729SShri Abhyankar for (j=0; j<nz; j++) { 27417542729SShri Abhyankar x = rtmp+9*ajtmp[j]; 27517542729SShri Abhyankar x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = 0.0; 27617542729SShri Abhyankar } 27717542729SShri Abhyankar /* load in initial (unfactored row) */ 27817542729SShri Abhyankar nz = ai[i+1] - ai[i]; 27917542729SShri Abhyankar ajtmpold = aj + ai[i]; 28017542729SShri Abhyankar v = aa + 9*ai[i]; 28117542729SShri Abhyankar for (j=0; j<nz; j++) { 28217542729SShri Abhyankar x = rtmp+9*ajtmpold[j]; 28317542729SShri Abhyankar x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 28417542729SShri Abhyankar x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8]; 28517542729SShri Abhyankar v += 9; 28617542729SShri Abhyankar } 28717542729SShri Abhyankar row = *ajtmp++; 28817542729SShri Abhyankar while (row < i) { 28917542729SShri Abhyankar pc = rtmp + 9*row; 29017542729SShri Abhyankar p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 29117542729SShri Abhyankar p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8]; 29217542729SShri Abhyankar if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || 29317542729SShri Abhyankar p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) { 29417542729SShri Abhyankar pv = ba + 9*diag_offset[row]; 29517542729SShri Abhyankar pj = bj + diag_offset[row] + 1; 29617542729SShri Abhyankar x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 29717542729SShri Abhyankar x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 29817542729SShri Abhyankar pc[0] = m1 = p1*x1 + p4*x2 + p7*x3; 29917542729SShri Abhyankar pc[1] = m2 = p2*x1 + p5*x2 + p8*x3; 30017542729SShri Abhyankar pc[2] = m3 = p3*x1 + p6*x2 + p9*x3; 30117542729SShri Abhyankar 30217542729SShri Abhyankar pc[3] = m4 = p1*x4 + p4*x5 + p7*x6; 30317542729SShri Abhyankar pc[4] = m5 = p2*x4 + p5*x5 + p8*x6; 30417542729SShri Abhyankar pc[5] = m6 = p3*x4 + p6*x5 + p9*x6; 30517542729SShri Abhyankar 30617542729SShri Abhyankar pc[6] = m7 = p1*x7 + p4*x8 + p7*x9; 30717542729SShri Abhyankar pc[7] = m8 = p2*x7 + p5*x8 + p8*x9; 30817542729SShri Abhyankar pc[8] = m9 = p3*x7 + p6*x8 + p9*x9; 30917542729SShri Abhyankar 31017542729SShri Abhyankar nz = bi[row+1] - diag_offset[row] - 1; 31117542729SShri Abhyankar pv += 9; 31217542729SShri Abhyankar for (j=0; j<nz; j++) { 31317542729SShri Abhyankar x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 31417542729SShri Abhyankar x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 31517542729SShri Abhyankar x = rtmp + 9*pj[j]; 31617542729SShri Abhyankar x[0] -= m1*x1 + m4*x2 + m7*x3; 31717542729SShri Abhyankar x[1] -= m2*x1 + m5*x2 + m8*x3; 31817542729SShri Abhyankar x[2] -= m3*x1 + m6*x2 + m9*x3; 31917542729SShri Abhyankar 32017542729SShri Abhyankar x[3] -= m1*x4 + m4*x5 + m7*x6; 32117542729SShri Abhyankar x[4] -= m2*x4 + m5*x5 + m8*x6; 32217542729SShri Abhyankar x[5] -= m3*x4 + m6*x5 + m9*x6; 32317542729SShri Abhyankar 32417542729SShri Abhyankar x[6] -= m1*x7 + m4*x8 + m7*x9; 32517542729SShri Abhyankar x[7] -= m2*x7 + m5*x8 + m8*x9; 32617542729SShri Abhyankar x[8] -= m3*x7 + m6*x8 + m9*x9; 32717542729SShri Abhyankar pv += 9; 32817542729SShri Abhyankar } 32917542729SShri Abhyankar ierr = PetscLogFlops(54.0*nz+36.0);CHKERRQ(ierr); 33017542729SShri Abhyankar } 33117542729SShri Abhyankar row = *ajtmp++; 33217542729SShri Abhyankar } 33317542729SShri Abhyankar /* finished row so stick it into b->a */ 33417542729SShri Abhyankar pv = ba + 9*bi[i]; 33517542729SShri Abhyankar pj = bj + bi[i]; 33617542729SShri Abhyankar nz = bi[i+1] - bi[i]; 33717542729SShri Abhyankar for (j=0; j<nz; j++) { 33817542729SShri Abhyankar x = rtmp+9*pj[j]; 33917542729SShri Abhyankar pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 34017542729SShri Abhyankar pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8]; 34117542729SShri Abhyankar pv += 9; 34217542729SShri Abhyankar } 34317542729SShri Abhyankar /* invert diagonal block */ 34417542729SShri Abhyankar w = ba + 9*diag_offset[i]; 345*a455e926SHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 346*a455e926SHong Zhang ierr = PetscKernel_A_gets_inverse_A_3(w,shift,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr); 347922032d7SHong Zhang if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 34817542729SShri Abhyankar } 34917542729SShri Abhyankar 35017542729SShri Abhyankar ierr = PetscFree(rtmp);CHKERRQ(ierr); 35126fbe8dcSKarl Rupp 35206e38f1dSHong Zhang C->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering_inplace; 35306e38f1dSHong Zhang C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering_inplace; 35417542729SShri Abhyankar C->assembled = PETSC_TRUE; 35526fbe8dcSKarl Rupp 356766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*3*3*3*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 35717542729SShri Abhyankar PetscFunctionReturn(0); 35817542729SShri Abhyankar } 35917542729SShri Abhyankar 36017542729SShri Abhyankar /* 3614dd39f65SShri Abhyankar MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering - 3624dd39f65SShri Abhyankar copied from MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace() 36317542729SShri Abhyankar */ 36417542729SShri Abhyankar #undef __FUNCT__ 3654dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering" 3664dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info) 36717542729SShri Abhyankar { 36817542729SShri Abhyankar Mat C =B; 36917542729SShri Abhyankar Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data; 37017542729SShri Abhyankar PetscErrorCode ierr; 371bbd65245SShri Abhyankar PetscInt i,j,k,nz,nzL,row; 372bbd65245SShri Abhyankar const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 373bbd65245SShri Abhyankar const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2; 37417542729SShri Abhyankar MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 375bbd65245SShri Abhyankar PetscInt flg; 376182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 377*a455e926SHong Zhang PetscBool allowzeropivot,zeropivotdetected; 37817542729SShri Abhyankar 37917542729SShri Abhyankar PetscFunctionBegin; 38017542729SShri Abhyankar /* generate work space needed by the factorization */ 381dcca6d9dSJed Brown ierr = PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);CHKERRQ(ierr); 38217542729SShri Abhyankar ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 38317542729SShri Abhyankar 38417542729SShri Abhyankar for (i=0; i<n; i++) { 38517542729SShri Abhyankar /* zero rtmp */ 38617542729SShri Abhyankar /* L part */ 38717542729SShri Abhyankar nz = bi[i+1] - bi[i]; 38817542729SShri Abhyankar bjtmp = bj + bi[i]; 38917542729SShri Abhyankar for (j=0; j<nz; j++) { 39017542729SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 39117542729SShri Abhyankar } 39217542729SShri Abhyankar 39317542729SShri Abhyankar /* U part */ 394b2b2dd24SShri Abhyankar nz = bdiag[i] - bdiag[i+1]; 395b2b2dd24SShri Abhyankar bjtmp = bj + bdiag[i+1] + 1; 396b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 397b2b2dd24SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 398b2b2dd24SShri Abhyankar } 399b2b2dd24SShri Abhyankar 400b2b2dd24SShri Abhyankar /* load in initial (unfactored row) */ 401b2b2dd24SShri Abhyankar nz = ai[i+1] - ai[i]; 402b2b2dd24SShri Abhyankar ajtmp = aj + ai[i]; 403b2b2dd24SShri Abhyankar v = aa + bs2*ai[i]; 404b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 405b2b2dd24SShri Abhyankar ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 406b2b2dd24SShri Abhyankar } 407b2b2dd24SShri Abhyankar 408b2b2dd24SShri Abhyankar /* elimination */ 409b2b2dd24SShri Abhyankar bjtmp = bj + bi[i]; 410b2b2dd24SShri Abhyankar nzL = bi[i+1] - bi[i]; 411b2b2dd24SShri Abhyankar for (k=0; k<nzL; k++) { 412b2b2dd24SShri Abhyankar row = bjtmp[k]; 413b2b2dd24SShri Abhyankar pc = rtmp + bs2*row; 414c35f09e5SBarry Smith for (flg=0,j=0; j<bs2; j++) { 415c35f09e5SBarry Smith if (pc[j]!=0.0) { 416c35f09e5SBarry Smith flg = 1; 417c35f09e5SBarry Smith break; 418c35f09e5SBarry Smith } 419c35f09e5SBarry Smith } 420b2b2dd24SShri Abhyankar if (flg) { 421b2b2dd24SShri Abhyankar pv = b->a + bs2*bdiag[row]; 42296b95a6bSBarry Smith /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 42396b95a6bSBarry Smith ierr = PetscKernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr); 424b2b2dd24SShri Abhyankar 425b2b2dd24SShri Abhyankar pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */ 426b2b2dd24SShri Abhyankar pv = b->a + bs2*(bdiag[row+1]+1); 427b2b2dd24SShri Abhyankar nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */ 428b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 42996b95a6bSBarry Smith /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 430b2b2dd24SShri Abhyankar /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 431b2b2dd24SShri Abhyankar v = rtmp + bs2*pj[j]; 43296b95a6bSBarry Smith ierr = PetscKernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr); 433b2b2dd24SShri Abhyankar pv += bs2; 434b2b2dd24SShri Abhyankar } 435b2b2dd24SShri Abhyankar ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 436b2b2dd24SShri Abhyankar } 437b2b2dd24SShri Abhyankar } 438b2b2dd24SShri Abhyankar 439b2b2dd24SShri Abhyankar /* finished row so stick it into b->a */ 440b2b2dd24SShri Abhyankar /* L part */ 441b2b2dd24SShri Abhyankar pv = b->a + bs2*bi[i]; 442b2b2dd24SShri Abhyankar pj = b->j + bi[i]; 443b2b2dd24SShri Abhyankar nz = bi[i+1] - bi[i]; 444b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 445b2b2dd24SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 446b2b2dd24SShri Abhyankar } 447b2b2dd24SShri Abhyankar 448b2b2dd24SShri Abhyankar /* Mark diagonal and invert diagonal for simplier triangular solves */ 449b2b2dd24SShri Abhyankar pv = b->a + bs2*bdiag[i]; 450b2b2dd24SShri Abhyankar pj = b->j + bdiag[i]; 451b2b2dd24SShri Abhyankar ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 45296b95a6bSBarry Smith /* ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 453*a455e926SHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 454*a455e926SHong Zhang ierr = PetscKernel_A_gets_inverse_A_3(pv,shift,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr); 455922032d7SHong Zhang if (zeropivotdetected) B->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 456b2b2dd24SShri Abhyankar 457b2b2dd24SShri Abhyankar /* U part */ 458b2b2dd24SShri Abhyankar pv = b->a + bs2*(bdiag[i+1]+1); 459b2b2dd24SShri Abhyankar pj = b->j + bdiag[i+1]+1; 460b2b2dd24SShri Abhyankar nz = bdiag[i] - bdiag[i+1] - 1; 461b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 462b2b2dd24SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 463b2b2dd24SShri Abhyankar } 464b2b2dd24SShri Abhyankar } 465fca92195SBarry Smith ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); 46626fbe8dcSKarl Rupp 4674dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering; 4689f5c0bcdSBarry Smith C->ops->forwardsolve = MatForwardSolve_SeqBAIJ_3_NaturalOrdering; 4699f5c0bcdSBarry Smith C->ops->backwardsolve = MatBackwardSolve_SeqBAIJ_3_NaturalOrdering; 4704dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering; 471b2b2dd24SShri Abhyankar C->assembled = PETSC_TRUE; 47226fbe8dcSKarl Rupp 473766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*3*3*3*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 474b2b2dd24SShri Abhyankar PetscFunctionReturn(0); 475b2b2dd24SShri Abhyankar } 476b2b2dd24SShri Abhyankar 477