183287d42SBarry Smith /* 283287d42SBarry Smith Factorization code for BAIJ format. 383287d42SBarry Smith */ 483287d42SBarry Smith #include "src/mat/impls/baij/seq/baij.h" 583287d42SBarry Smith #include "src/inline/ilu.h" 683287d42SBarry Smith 783287d42SBarry Smith /* ------------------------------------------------------------*/ 883287d42SBarry Smith /* 983287d42SBarry Smith Version for when blocks are 4 by 4 1083287d42SBarry Smith */ 114a2ae208SSatish Balay #undef __FUNCT__ 124a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_4" 13dfbe8321SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4(Mat A,Mat *B) 1483287d42SBarry Smith { 1583287d42SBarry Smith Mat C = *B; 1683287d42SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 1783287d42SBarry Smith IS isrow = b->row,isicol = b->icol; 18*6849ba73SBarry Smith PetscErrorCode ierr; 19*6849ba73SBarry Smith int *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j; 2083287d42SBarry Smith int *ajtmpold,*ajtmp,nz,row; 2183287d42SBarry Smith int *diag_offset = b->diag,idx,*ai=a->i,*aj=a->j,*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,x10,x11,x12,x13,x14,x15,x16; 2583287d42SBarry Smith MatScalar p10,p11,p12,p13,p14,p15,p16,m10,m11,m12; 2683287d42SBarry Smith MatScalar m13,m14,m15,m16; 2783287d42SBarry Smith MatScalar *ba = b->a,*aa = a->a; 28bcd9e38bSBarry Smith PetscTruth pivotinblocks = b->pivotinblocks; 2983287d42SBarry Smith 3083287d42SBarry Smith PetscFunctionBegin; 3183287d42SBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 3283287d42SBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 33b0a32e0cSBarry Smith ierr = PetscMalloc(16*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 3483287d42SBarry Smith 3583287d42SBarry Smith for (i=0; i<n; i++) { 3683287d42SBarry Smith nz = bi[i+1] - bi[i]; 3783287d42SBarry Smith ajtmp = bj + bi[i]; 3883287d42SBarry Smith for (j=0; j<nz; j++) { 3983287d42SBarry Smith x = rtmp+16*ajtmp[j]; 4083287d42SBarry Smith x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0; 4183287d42SBarry Smith x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = 0.0; 4283287d42SBarry Smith } 4383287d42SBarry Smith /* load in initial (unfactored row) */ 4483287d42SBarry Smith idx = r[i]; 4583287d42SBarry Smith nz = ai[idx+1] - ai[idx]; 4683287d42SBarry Smith ajtmpold = aj + ai[idx]; 4783287d42SBarry Smith v = aa + 16*ai[idx]; 4883287d42SBarry Smith for (j=0; j<nz; j++) { 4983287d42SBarry Smith x = rtmp+16*ic[ajtmpold[j]]; 5083287d42SBarry Smith x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 5183287d42SBarry Smith x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8]; 5283287d42SBarry Smith x[9] = v[9]; x[10] = v[10]; x[11] = v[11]; x[12] = v[12]; x[13] = v[13]; 5383287d42SBarry Smith x[14] = v[14]; x[15] = v[15]; 5483287d42SBarry Smith v += 16; 5583287d42SBarry Smith } 5683287d42SBarry Smith row = *ajtmp++; 5783287d42SBarry Smith while (row < i) { 5883287d42SBarry Smith pc = rtmp + 16*row; 5983287d42SBarry Smith p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 6083287d42SBarry Smith p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8]; 6183287d42SBarry Smith p10 = pc[9]; p11 = pc[10]; p12 = pc[11]; p13 = pc[12]; p14 = pc[13]; 6283287d42SBarry Smith p15 = pc[14]; p16 = pc[15]; 6383287d42SBarry Smith if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || 6483287d42SBarry Smith p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 || 6583287d42SBarry Smith p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0 6683287d42SBarry Smith || p16 != 0.0) { 6783287d42SBarry Smith pv = ba + 16*diag_offset[row]; 6883287d42SBarry Smith pj = bj + diag_offset[row] + 1; 6983287d42SBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 7083287d42SBarry Smith x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 7183287d42SBarry Smith x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; x14 = pv[13]; 7283287d42SBarry Smith x15 = pv[14]; x16 = pv[15]; 7383287d42SBarry Smith pc[0] = m1 = p1*x1 + p5*x2 + p9*x3 + p13*x4; 7483287d42SBarry Smith pc[1] = m2 = p2*x1 + p6*x2 + p10*x3 + p14*x4; 7583287d42SBarry Smith pc[2] = m3 = p3*x1 + p7*x2 + p11*x3 + p15*x4; 7683287d42SBarry Smith pc[3] = m4 = p4*x1 + p8*x2 + p12*x3 + p16*x4; 7783287d42SBarry Smith 7883287d42SBarry Smith pc[4] = m5 = p1*x5 + p5*x6 + p9*x7 + p13*x8; 7983287d42SBarry Smith pc[5] = m6 = p2*x5 + p6*x6 + p10*x7 + p14*x8; 8083287d42SBarry Smith pc[6] = m7 = p3*x5 + p7*x6 + p11*x7 + p15*x8; 8183287d42SBarry Smith pc[7] = m8 = p4*x5 + p8*x6 + p12*x7 + p16*x8; 8283287d42SBarry Smith 8383287d42SBarry Smith pc[8] = m9 = p1*x9 + p5*x10 + p9*x11 + p13*x12; 8483287d42SBarry Smith pc[9] = m10 = p2*x9 + p6*x10 + p10*x11 + p14*x12; 8583287d42SBarry Smith pc[10] = m11 = p3*x9 + p7*x10 + p11*x11 + p15*x12; 8683287d42SBarry Smith pc[11] = m12 = p4*x9 + p8*x10 + p12*x11 + p16*x12; 8783287d42SBarry Smith 8883287d42SBarry Smith pc[12] = m13 = p1*x13 + p5*x14 + p9*x15 + p13*x16; 8983287d42SBarry Smith pc[13] = m14 = p2*x13 + p6*x14 + p10*x15 + p14*x16; 9083287d42SBarry Smith pc[14] = m15 = p3*x13 + p7*x14 + p11*x15 + p15*x16; 9183287d42SBarry Smith pc[15] = m16 = p4*x13 + p8*x14 + p12*x15 + p16*x16; 9283287d42SBarry Smith 9383287d42SBarry Smith nz = bi[row+1] - diag_offset[row] - 1; 9483287d42SBarry Smith pv += 16; 9583287d42SBarry Smith for (j=0; j<nz; j++) { 9683287d42SBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 9783287d42SBarry Smith x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 9883287d42SBarry Smith x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; 9983287d42SBarry Smith x14 = pv[13]; x15 = pv[14]; x16 = pv[15]; 10083287d42SBarry Smith x = rtmp + 16*pj[j]; 10183287d42SBarry Smith x[0] -= m1*x1 + m5*x2 + m9*x3 + m13*x4; 10283287d42SBarry Smith x[1] -= m2*x1 + m6*x2 + m10*x3 + m14*x4; 10383287d42SBarry Smith x[2] -= m3*x1 + m7*x2 + m11*x3 + m15*x4; 10483287d42SBarry Smith x[3] -= m4*x1 + m8*x2 + m12*x3 + m16*x4; 10583287d42SBarry Smith 10683287d42SBarry Smith x[4] -= m1*x5 + m5*x6 + m9*x7 + m13*x8; 10783287d42SBarry Smith x[5] -= m2*x5 + m6*x6 + m10*x7 + m14*x8; 10883287d42SBarry Smith x[6] -= m3*x5 + m7*x6 + m11*x7 + m15*x8; 10983287d42SBarry Smith x[7] -= m4*x5 + m8*x6 + m12*x7 + m16*x8; 11083287d42SBarry Smith 11183287d42SBarry Smith x[8] -= m1*x9 + m5*x10 + m9*x11 + m13*x12; 11283287d42SBarry Smith x[9] -= m2*x9 + m6*x10 + m10*x11 + m14*x12; 11383287d42SBarry Smith x[10] -= m3*x9 + m7*x10 + m11*x11 + m15*x12; 11483287d42SBarry Smith x[11] -= m4*x9 + m8*x10 + m12*x11 + m16*x12; 11583287d42SBarry Smith 11683287d42SBarry Smith x[12] -= m1*x13 + m5*x14 + m9*x15 + m13*x16; 11783287d42SBarry Smith x[13] -= m2*x13 + m6*x14 + m10*x15 + m14*x16; 11883287d42SBarry Smith x[14] -= m3*x13 + m7*x14 + m11*x15 + m15*x16; 11983287d42SBarry Smith x[15] -= m4*x13 + m8*x14 + m12*x15 + m16*x16; 12083287d42SBarry Smith 12183287d42SBarry Smith pv += 16; 12283287d42SBarry Smith } 123b0a32e0cSBarry Smith PetscLogFlops(128*nz+112); 12483287d42SBarry Smith } 12583287d42SBarry Smith row = *ajtmp++; 12683287d42SBarry Smith } 12783287d42SBarry Smith /* finished row so stick it into b->a */ 12883287d42SBarry Smith pv = ba + 16*bi[i]; 12983287d42SBarry Smith pj = bj + bi[i]; 13083287d42SBarry Smith nz = bi[i+1] - bi[i]; 13183287d42SBarry Smith for (j=0; j<nz; j++) { 13283287d42SBarry Smith x = rtmp+16*pj[j]; 13383287d42SBarry Smith pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 13483287d42SBarry Smith pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8]; 13583287d42SBarry Smith pv[9] = x[9]; pv[10] = x[10]; pv[11] = x[11]; pv[12] = x[12]; 13683287d42SBarry Smith pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15]; 13783287d42SBarry Smith pv += 16; 13883287d42SBarry Smith } 13983287d42SBarry Smith /* invert diagonal block */ 14083287d42SBarry Smith w = ba + 16*diag_offset[i]; 141bcd9e38bSBarry Smith if (pivotinblocks) { 14283287d42SBarry Smith ierr = Kernel_A_gets_inverse_A_4(w);CHKERRQ(ierr); 143bcd9e38bSBarry Smith } else { 144bcd9e38bSBarry Smith ierr = Kernel_A_gets_inverse_A_4_nopivot(w);CHKERRQ(ierr); 145bcd9e38bSBarry Smith } 14683287d42SBarry Smith } 14783287d42SBarry Smith 14883287d42SBarry Smith ierr = PetscFree(rtmp);CHKERRQ(ierr); 14983287d42SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 15083287d42SBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 15183287d42SBarry Smith C->factor = FACTOR_LU; 15283287d42SBarry Smith C->assembled = PETSC_TRUE; 153b0a32e0cSBarry Smith PetscLogFlops(1.3333*64*b->mbs); /* from inverting diagonal blocks */ 15483287d42SBarry Smith PetscFunctionReturn(0); 15583287d42SBarry Smith } 156