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__ 134a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3" 140481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3(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); 115db4efbfdSBarry Smith C->ops->solve = MatSolve_SeqBAIJ_3; 116db4efbfdSBarry Smith C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3; 11783287d42SBarry Smith C->assembled = PETSC_TRUE; 118efee365bSSatish Balay ierr = PetscLogFlops(1.3333*27*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 11983287d42SBarry Smith PetscFunctionReturn(0); 12083287d42SBarry Smith } 121*17542729SShri Abhyankar 122*17542729SShri Abhyankar /* MatLUFactorNumeric_SeqBAIJ_3_newdatastruct - 123*17542729SShri Abhyankar copied from MatLUFactorNumeric_SeqBAIJ_N_newdatastruct() and manually re-implemented 124*17542729SShri Abhyankar Kernel_A_gets_A_times_B() 125*17542729SShri Abhyankar Kernel_A_gets_A_minus_B_times_C() 126*17542729SShri Abhyankar Kernel_A_gets_inverse_A() 127*17542729SShri Abhyankar */ 128*17542729SShri Abhyankar #undef __FUNCT__ 129*17542729SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_newdatastruct" 130*17542729SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_newdatastruct(Mat B,Mat A,const MatFactorInfo *info) 131*17542729SShri Abhyankar { 132*17542729SShri Abhyankar Mat C=B; 133*17542729SShri Abhyankar Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; 134*17542729SShri Abhyankar IS isrow = b->row,isicol = b->icol; 135*17542729SShri Abhyankar PetscErrorCode ierr; 136*17542729SShri Abhyankar const PetscInt *r,*ic,*ics; 137*17542729SShri Abhyankar PetscInt i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 138*17542729SShri Abhyankar PetscInt *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj; 139*17542729SShri Abhyankar MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 140*17542729SShri Abhyankar PetscInt bs2 = a->bs2,flg; 141*17542729SShri Abhyankar PetscReal shift = info->shiftinblocks; 142*17542729SShri Abhyankar 143*17542729SShri Abhyankar PetscFunctionBegin; 144*17542729SShri Abhyankar ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 145*17542729SShri Abhyankar ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 146*17542729SShri Abhyankar 147*17542729SShri Abhyankar /* generate work space needed by the factorization */ 148*17542729SShri Abhyankar ierr = PetscMalloc((bs2*n+bs2+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 149*17542729SShri Abhyankar mwork = rtmp + bs2*n; 150*17542729SShri Abhyankar ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 151*17542729SShri Abhyankar ics = ic; 152*17542729SShri Abhyankar 153*17542729SShri Abhyankar for (i=0; i<n; i++){ 154*17542729SShri Abhyankar /* zero rtmp */ 155*17542729SShri Abhyankar /* L part */ 156*17542729SShri Abhyankar nz = bi[i+1] - bi[i]; 157*17542729SShri Abhyankar bjtmp = bj + bi[i]; 158*17542729SShri Abhyankar for (j=0; j<nz; j++){ 159*17542729SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 160*17542729SShri Abhyankar } 161*17542729SShri Abhyankar 162*17542729SShri Abhyankar /* U part */ 163*17542729SShri Abhyankar nz = bi[2*n-i+1] - bi[2*n-i]; 164*17542729SShri Abhyankar bjtmp = bj + bi[2*n-i]; 165*17542729SShri Abhyankar for (j=0; j<nz; j++){ 166*17542729SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 167*17542729SShri Abhyankar } 168*17542729SShri Abhyankar 169*17542729SShri Abhyankar /* load in initial (unfactored row) */ 170*17542729SShri Abhyankar nz = ai[r[i]+1] - ai[r[i]]; 171*17542729SShri Abhyankar ajtmp = aj + ai[r[i]]; 172*17542729SShri Abhyankar v = aa + bs2*ai[r[i]]; 173*17542729SShri Abhyankar for (j=0; j<nz; j++) { 174*17542729SShri Abhyankar ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 175*17542729SShri Abhyankar } 176*17542729SShri Abhyankar 177*17542729SShri Abhyankar /* elimination */ 178*17542729SShri Abhyankar bjtmp = bj + bi[i]; 179*17542729SShri Abhyankar row = *bjtmp++; 180*17542729SShri Abhyankar nzL = bi[i+1] - bi[i]; 181*17542729SShri Abhyankar k = 0; 182*17542729SShri Abhyankar while (k < nzL) { 183*17542729SShri Abhyankar pc = rtmp + bs2*row; 184*17542729SShri Abhyankar for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }} 185*17542729SShri Abhyankar if (flg) { 186*17542729SShri Abhyankar pv = b->a + bs2*bdiag[row]; 187*17542729SShri Abhyankar /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 188*17542729SShri Abhyankar ierr = Kernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr); 189*17542729SShri Abhyankar 190*17542729SShri Abhyankar pj = b->j + bi[2*n-row]; /* begining of U(row,:) */ 191*17542729SShri Abhyankar pv = b->a + bs2*bi[2*n-row]; 192*17542729SShri Abhyankar nz = bi[2*n-row+1] - bi[2*n-row] - 1; /* num of entries inU(row,:), excluding diag */ 193*17542729SShri Abhyankar for (j=0; j<nz; j++) { 194*17542729SShri Abhyankar /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 195*17542729SShri Abhyankar /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 196*17542729SShri Abhyankar v = rtmp + bs2*pj[j]; 197*17542729SShri Abhyankar ierr = Kernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr); 198*17542729SShri Abhyankar pv += bs2; 199*17542729SShri Abhyankar } 200*17542729SShri Abhyankar ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 201*17542729SShri Abhyankar } 202*17542729SShri Abhyankar row = *bjtmp++; k++; 203*17542729SShri Abhyankar } 204*17542729SShri Abhyankar 205*17542729SShri Abhyankar /* finished row so stick it into b->a */ 206*17542729SShri Abhyankar /* L part */ 207*17542729SShri Abhyankar pv = b->a + bs2*bi[i] ; 208*17542729SShri Abhyankar pj = b->j + bi[i] ; 209*17542729SShri Abhyankar nz = bi[i+1] - bi[i]; 210*17542729SShri Abhyankar for (j=0; j<nz; j++) { 211*17542729SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 212*17542729SShri Abhyankar } 213*17542729SShri Abhyankar 214*17542729SShri Abhyankar /* Mark diagonal and invert diagonal for simplier triangular solves */ 215*17542729SShri Abhyankar pv = b->a + bs2*bdiag[i]; 216*17542729SShri Abhyankar pj = b->j + bdiag[i]; 217*17542729SShri Abhyankar ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 218*17542729SShri Abhyankar /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 219*17542729SShri Abhyankar ierr = Kernel_A_gets_inverse_A_3(pv,shift);CHKERRQ(ierr); 220*17542729SShri Abhyankar 221*17542729SShri Abhyankar /* U part */ 222*17542729SShri Abhyankar pv = b->a + bs2*bi[2*n-i]; 223*17542729SShri Abhyankar pj = b->j + bi[2*n-i]; 224*17542729SShri Abhyankar nz = bi[2*n-i+1] - bi[2*n-i] - 1; 225*17542729SShri Abhyankar for (j=0; j<nz; j++){ 226*17542729SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 227*17542729SShri Abhyankar } 228*17542729SShri Abhyankar } 229*17542729SShri Abhyankar 230*17542729SShri Abhyankar ierr = PetscFree(rtmp);CHKERRQ(ierr); 231*17542729SShri Abhyankar ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 232*17542729SShri Abhyankar ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 233*17542729SShri Abhyankar 234*17542729SShri Abhyankar C->assembled = PETSC_TRUE; 235*17542729SShri Abhyankar ierr = PetscLogFlops(1.3333*bs2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 236*17542729SShri Abhyankar PetscFunctionReturn(0); 237*17542729SShri Abhyankar } 238*17542729SShri Abhyankar 239*17542729SShri Abhyankar #undef __FUNCT__ 240*17542729SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering" 241*17542729SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info) 242*17542729SShri Abhyankar { 243*17542729SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 244*17542729SShri Abhyankar PetscErrorCode ierr; 245*17542729SShri Abhyankar PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 246*17542729SShri Abhyankar PetscInt *ajtmpold,*ajtmp,nz,row; 247*17542729SShri Abhyankar PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; 248*17542729SShri Abhyankar MatScalar *pv,*v,*rtmp,*pc,*w,*x; 249*17542729SShri Abhyankar MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4; 250*17542729SShri Abhyankar MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9; 251*17542729SShri Abhyankar MatScalar *ba = b->a,*aa = a->a; 252*17542729SShri Abhyankar PetscReal shift = info->shiftinblocks; 253*17542729SShri Abhyankar 254*17542729SShri Abhyankar PetscFunctionBegin; 255*17542729SShri Abhyankar ierr = PetscMalloc(9*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 256*17542729SShri Abhyankar 257*17542729SShri Abhyankar for (i=0; i<n; i++) { 258*17542729SShri Abhyankar nz = bi[i+1] - bi[i]; 259*17542729SShri Abhyankar ajtmp = bj + bi[i]; 260*17542729SShri Abhyankar for (j=0; j<nz; j++) { 261*17542729SShri Abhyankar x = rtmp+9*ajtmp[j]; 262*17542729SShri Abhyankar x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = 0.0; 263*17542729SShri Abhyankar } 264*17542729SShri Abhyankar /* load in initial (unfactored row) */ 265*17542729SShri Abhyankar nz = ai[i+1] - ai[i]; 266*17542729SShri Abhyankar ajtmpold = aj + ai[i]; 267*17542729SShri Abhyankar v = aa + 9*ai[i]; 268*17542729SShri Abhyankar for (j=0; j<nz; j++) { 269*17542729SShri Abhyankar x = rtmp+9*ajtmpold[j]; 270*17542729SShri Abhyankar x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 271*17542729SShri Abhyankar x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8]; 272*17542729SShri Abhyankar v += 9; 273*17542729SShri Abhyankar } 274*17542729SShri Abhyankar row = *ajtmp++; 275*17542729SShri Abhyankar while (row < i) { 276*17542729SShri Abhyankar pc = rtmp + 9*row; 277*17542729SShri Abhyankar p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 278*17542729SShri Abhyankar p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8]; 279*17542729SShri Abhyankar if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || 280*17542729SShri Abhyankar p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) { 281*17542729SShri Abhyankar pv = ba + 9*diag_offset[row]; 282*17542729SShri Abhyankar pj = bj + diag_offset[row] + 1; 283*17542729SShri Abhyankar x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 284*17542729SShri Abhyankar x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 285*17542729SShri Abhyankar pc[0] = m1 = p1*x1 + p4*x2 + p7*x3; 286*17542729SShri Abhyankar pc[1] = m2 = p2*x1 + p5*x2 + p8*x3; 287*17542729SShri Abhyankar pc[2] = m3 = p3*x1 + p6*x2 + p9*x3; 288*17542729SShri Abhyankar 289*17542729SShri Abhyankar pc[3] = m4 = p1*x4 + p4*x5 + p7*x6; 290*17542729SShri Abhyankar pc[4] = m5 = p2*x4 + p5*x5 + p8*x6; 291*17542729SShri Abhyankar pc[5] = m6 = p3*x4 + p6*x5 + p9*x6; 292*17542729SShri Abhyankar 293*17542729SShri Abhyankar pc[6] = m7 = p1*x7 + p4*x8 + p7*x9; 294*17542729SShri Abhyankar pc[7] = m8 = p2*x7 + p5*x8 + p8*x9; 295*17542729SShri Abhyankar pc[8] = m9 = p3*x7 + p6*x8 + p9*x9; 296*17542729SShri Abhyankar 297*17542729SShri Abhyankar nz = bi[row+1] - diag_offset[row] - 1; 298*17542729SShri Abhyankar pv += 9; 299*17542729SShri Abhyankar for (j=0; j<nz; j++) { 300*17542729SShri Abhyankar x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 301*17542729SShri Abhyankar x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 302*17542729SShri Abhyankar x = rtmp + 9*pj[j]; 303*17542729SShri Abhyankar x[0] -= m1*x1 + m4*x2 + m7*x3; 304*17542729SShri Abhyankar x[1] -= m2*x1 + m5*x2 + m8*x3; 305*17542729SShri Abhyankar x[2] -= m3*x1 + m6*x2 + m9*x3; 306*17542729SShri Abhyankar 307*17542729SShri Abhyankar x[3] -= m1*x4 + m4*x5 + m7*x6; 308*17542729SShri Abhyankar x[4] -= m2*x4 + m5*x5 + m8*x6; 309*17542729SShri Abhyankar x[5] -= m3*x4 + m6*x5 + m9*x6; 310*17542729SShri Abhyankar 311*17542729SShri Abhyankar x[6] -= m1*x7 + m4*x8 + m7*x9; 312*17542729SShri Abhyankar x[7] -= m2*x7 + m5*x8 + m8*x9; 313*17542729SShri Abhyankar x[8] -= m3*x7 + m6*x8 + m9*x9; 314*17542729SShri Abhyankar pv += 9; 315*17542729SShri Abhyankar } 316*17542729SShri Abhyankar ierr = PetscLogFlops(54.0*nz+36.0);CHKERRQ(ierr); 317*17542729SShri Abhyankar } 318*17542729SShri Abhyankar row = *ajtmp++; 319*17542729SShri Abhyankar } 320*17542729SShri Abhyankar /* finished row so stick it into b->a */ 321*17542729SShri Abhyankar pv = ba + 9*bi[i]; 322*17542729SShri Abhyankar pj = bj + bi[i]; 323*17542729SShri Abhyankar nz = bi[i+1] - bi[i]; 324*17542729SShri Abhyankar for (j=0; j<nz; j++) { 325*17542729SShri Abhyankar x = rtmp+9*pj[j]; 326*17542729SShri Abhyankar pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 327*17542729SShri Abhyankar pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8]; 328*17542729SShri Abhyankar pv += 9; 329*17542729SShri Abhyankar } 330*17542729SShri Abhyankar /* invert diagonal block */ 331*17542729SShri Abhyankar w = ba + 9*diag_offset[i]; 332*17542729SShri Abhyankar ierr = Kernel_A_gets_inverse_A_3(w,shift);CHKERRQ(ierr); 333*17542729SShri Abhyankar } 334*17542729SShri Abhyankar 335*17542729SShri Abhyankar ierr = PetscFree(rtmp);CHKERRQ(ierr); 336*17542729SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering; 337*17542729SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering; 338*17542729SShri Abhyankar C->assembled = PETSC_TRUE; 339*17542729SShri Abhyankar ierr = PetscLogFlops(1.3333*27*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 340*17542729SShri Abhyankar PetscFunctionReturn(0); 341*17542729SShri Abhyankar } 342*17542729SShri Abhyankar 343*17542729SShri Abhyankar /* 344*17542729SShri Abhyankar MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_newdatastruct - 345*17542729SShri Abhyankar copied from MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_newdatastruct() 346*17542729SShri Abhyankar */ 347*17542729SShri Abhyankar #undef __FUNCT__ 348*17542729SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_newdatastruct" 349*17542729SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_newdatastruct(Mat B,Mat A,const MatFactorInfo *info) 350*17542729SShri Abhyankar { 351*17542729SShri Abhyankar Mat C=B; 352*17542729SShri Abhyankar Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; 353*17542729SShri Abhyankar PetscErrorCode ierr; 354*17542729SShri Abhyankar PetscInt i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 355*17542729SShri Abhyankar PetscInt *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj; 356*17542729SShri Abhyankar MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 357*17542729SShri Abhyankar PetscInt bs2 = a->bs2,flg; 358*17542729SShri Abhyankar PetscReal shift = info->shiftinblocks; 359*17542729SShri Abhyankar 360*17542729SShri Abhyankar PetscFunctionBegin; 361*17542729SShri Abhyankar /* generate work space needed by the factorization */ 362*17542729SShri Abhyankar ierr = PetscMalloc((bs2*n+bs2+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 363*17542729SShri Abhyankar mwork = rtmp + bs2*n; 364*17542729SShri Abhyankar ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 365*17542729SShri Abhyankar 366*17542729SShri Abhyankar for (i=0; i<n; i++){ 367*17542729SShri Abhyankar /* zero rtmp */ 368*17542729SShri Abhyankar /* L part */ 369*17542729SShri Abhyankar nz = bi[i+1] - bi[i]; 370*17542729SShri Abhyankar bjtmp = bj + bi[i]; 371*17542729SShri Abhyankar for (j=0; j<nz; j++){ 372*17542729SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 373*17542729SShri Abhyankar } 374*17542729SShri Abhyankar 375*17542729SShri Abhyankar /* U part */ 376*17542729SShri Abhyankar nz = bi[2*n-i+1] - bi[2*n-i]; 377*17542729SShri Abhyankar bjtmp = bj + bi[2*n-i]; 378*17542729SShri Abhyankar for (j=0; j<nz; j++){ 379*17542729SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 380*17542729SShri Abhyankar } 381*17542729SShri Abhyankar 382*17542729SShri Abhyankar /* load in initial (unfactored row) */ 383*17542729SShri Abhyankar nz = ai[i+1] - ai[i]; 384*17542729SShri Abhyankar ajtmp = aj + ai[i]; 385*17542729SShri Abhyankar v = aa + bs2*ai[i]; 386*17542729SShri Abhyankar for (j=0; j<nz; j++) { 387*17542729SShri Abhyankar ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 388*17542729SShri Abhyankar } 389*17542729SShri Abhyankar 390*17542729SShri Abhyankar /* elimination */ 391*17542729SShri Abhyankar bjtmp = bj + bi[i]; 392*17542729SShri Abhyankar row = *bjtmp++; 393*17542729SShri Abhyankar nzL = bi[i+1] - bi[i]; 394*17542729SShri Abhyankar k = 0; 395*17542729SShri Abhyankar while (k < nzL) { 396*17542729SShri Abhyankar pc = rtmp + bs2*row; 397*17542729SShri Abhyankar for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }} 398*17542729SShri Abhyankar if (flg) { 399*17542729SShri Abhyankar pv = b->a + bs2*bdiag[row]; 400*17542729SShri Abhyankar /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 401*17542729SShri Abhyankar ierr = Kernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr); 402*17542729SShri Abhyankar 403*17542729SShri Abhyankar pj = b->j + bi[2*n-row]; /* begining of U(row,:) */ 404*17542729SShri Abhyankar pv = b->a + bs2*bi[2*n-row]; 405*17542729SShri Abhyankar nz = bi[2*n-row+1] - bi[2*n-row] - 1; /* num of entries inU(row,:), excluding diag */ 406*17542729SShri Abhyankar for (j=0; j<nz; j++) { 407*17542729SShri Abhyankar /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 408*17542729SShri Abhyankar /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 409*17542729SShri Abhyankar v = rtmp + bs2*pj[j]; 410*17542729SShri Abhyankar ierr = Kernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr); 411*17542729SShri Abhyankar pv += bs2; 412*17542729SShri Abhyankar } 413*17542729SShri Abhyankar ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 414*17542729SShri Abhyankar } 415*17542729SShri Abhyankar row = *bjtmp++; k++; 416*17542729SShri Abhyankar } 417*17542729SShri Abhyankar 418*17542729SShri Abhyankar /* finished row so stick it into b->a */ 419*17542729SShri Abhyankar /* L part */ 420*17542729SShri Abhyankar pv = b->a + bs2*bi[i] ; 421*17542729SShri Abhyankar pj = b->j + bi[i] ; 422*17542729SShri Abhyankar nz = bi[i+1] - bi[i]; 423*17542729SShri Abhyankar for (j=0; j<nz; j++) { 424*17542729SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 425*17542729SShri Abhyankar } 426*17542729SShri Abhyankar 427*17542729SShri Abhyankar /* Mark diagonal and invert diagonal for simplier triangular solves */ 428*17542729SShri Abhyankar pv = b->a + bs2*bdiag[i]; 429*17542729SShri Abhyankar pj = b->j + bdiag[i]; 430*17542729SShri Abhyankar ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 431*17542729SShri Abhyankar /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 432*17542729SShri Abhyankar ierr = Kernel_A_gets_inverse_A_3(pv,shift);CHKERRQ(ierr); 433*17542729SShri Abhyankar 434*17542729SShri Abhyankar /* U part */ 435*17542729SShri Abhyankar pv = b->a + bs2*bi[2*n-i]; 436*17542729SShri Abhyankar pj = b->j + bi[2*n-i]; 437*17542729SShri Abhyankar nz = bi[2*n-i+1] - bi[2*n-i] - 1; 438*17542729SShri Abhyankar for (j=0; j<nz; j++){ 439*17542729SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 440*17542729SShri Abhyankar } 441*17542729SShri Abhyankar } 442*17542729SShri Abhyankar 443*17542729SShri Abhyankar ierr = PetscFree(rtmp);CHKERRQ(ierr); 444*17542729SShri Abhyankar C->assembled = PETSC_TRUE; 445*17542729SShri Abhyankar ierr = PetscLogFlops(1.3333*bs2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 446*17542729SShri Abhyankar PetscFunctionReturn(0); 447*17542729SShri Abhyankar } 448