1be1d678aSKris Buschelman 283287d42SBarry Smith /* 383287d42SBarry Smith Factorization code for BAIJ format. 483287d42SBarry Smith */ 5c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h> 6c6db04a5SJed Brown #include <../src/mat/blockinvert.h> 783287d42SBarry Smith 883287d42SBarry Smith /* ------------------------------------------------------------*/ 983287d42SBarry Smith /* 1083287d42SBarry Smith Version for when blocks are 4 by 4 1183287d42SBarry Smith */ 124a2ae208SSatish Balay #undef __FUNCT__ 1306e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_4_inplace" 1406e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_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; 22690b6cddSBarry Smith PetscInt *diag_offset = b->diag,idx,*ai=a->i,*aj=a->j,*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,x10,x11,x12,x13,x14,x15,x16; 2683287d42SBarry Smith MatScalar p10,p11,p12,p13,p14,p15,p16,m10,m11,m12; 2783287d42SBarry Smith MatScalar m13,m14,m15,m16; 2883287d42SBarry Smith MatScalar *ba = b->a,*aa = a->a; 29ace3abfcSBarry Smith PetscBool pivotinblocks = b->pivotinblocks; 30182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 3183287d42SBarry Smith 3283287d42SBarry Smith PetscFunctionBegin; 3383287d42SBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 3483287d42SBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 35b0a32e0cSBarry Smith ierr = PetscMalloc(16*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 3683287d42SBarry Smith 3783287d42SBarry Smith for (i=0; i<n; i++) { 3883287d42SBarry Smith nz = bi[i+1] - bi[i]; 3983287d42SBarry Smith ajtmp = bj + bi[i]; 4083287d42SBarry Smith for (j=0; j<nz; j++) { 4183287d42SBarry Smith x = rtmp+16*ajtmp[j]; 4283287d42SBarry Smith x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0; 4383287d42SBarry Smith x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = 0.0; 4483287d42SBarry Smith } 4583287d42SBarry Smith /* load in initial (unfactored row) */ 4683287d42SBarry Smith idx = r[i]; 4783287d42SBarry Smith nz = ai[idx+1] - ai[idx]; 4883287d42SBarry Smith ajtmpold = aj + ai[idx]; 4983287d42SBarry Smith v = aa + 16*ai[idx]; 5083287d42SBarry Smith for (j=0; j<nz; j++) { 5183287d42SBarry Smith x = rtmp+16*ic[ajtmpold[j]]; 5283287d42SBarry Smith x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 5383287d42SBarry Smith x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8]; 5483287d42SBarry Smith x[9] = v[9]; x[10] = v[10]; x[11] = v[11]; x[12] = v[12]; x[13] = v[13]; 5583287d42SBarry Smith x[14] = v[14]; x[15] = v[15]; 5683287d42SBarry Smith v += 16; 5783287d42SBarry Smith } 5883287d42SBarry Smith row = *ajtmp++; 5983287d42SBarry Smith while (row < i) { 6083287d42SBarry Smith pc = rtmp + 16*row; 6183287d42SBarry Smith p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 6283287d42SBarry Smith p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8]; 6383287d42SBarry Smith p10 = pc[9]; p11 = pc[10]; p12 = pc[11]; p13 = pc[12]; p14 = pc[13]; 6483287d42SBarry Smith p15 = pc[14]; p16 = pc[15]; 6583287d42SBarry Smith if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || 6683287d42SBarry Smith p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 || 6783287d42SBarry Smith p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0 6883287d42SBarry Smith || p16 != 0.0) { 6983287d42SBarry Smith pv = ba + 16*diag_offset[row]; 7083287d42SBarry Smith pj = bj + diag_offset[row] + 1; 7183287d42SBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 7283287d42SBarry Smith x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 7383287d42SBarry Smith x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; x14 = pv[13]; 7483287d42SBarry Smith x15 = pv[14]; x16 = pv[15]; 7583287d42SBarry Smith pc[0] = m1 = p1*x1 + p5*x2 + p9*x3 + p13*x4; 7683287d42SBarry Smith pc[1] = m2 = p2*x1 + p6*x2 + p10*x3 + p14*x4; 7783287d42SBarry Smith pc[2] = m3 = p3*x1 + p7*x2 + p11*x3 + p15*x4; 7883287d42SBarry Smith pc[3] = m4 = p4*x1 + p8*x2 + p12*x3 + p16*x4; 7983287d42SBarry Smith 8083287d42SBarry Smith pc[4] = m5 = p1*x5 + p5*x6 + p9*x7 + p13*x8; 8183287d42SBarry Smith pc[5] = m6 = p2*x5 + p6*x6 + p10*x7 + p14*x8; 8283287d42SBarry Smith pc[6] = m7 = p3*x5 + p7*x6 + p11*x7 + p15*x8; 8383287d42SBarry Smith pc[7] = m8 = p4*x5 + p8*x6 + p12*x7 + p16*x8; 8483287d42SBarry Smith 8583287d42SBarry Smith pc[8] = m9 = p1*x9 + p5*x10 + p9*x11 + p13*x12; 8683287d42SBarry Smith pc[9] = m10 = p2*x9 + p6*x10 + p10*x11 + p14*x12; 8783287d42SBarry Smith pc[10] = m11 = p3*x9 + p7*x10 + p11*x11 + p15*x12; 8883287d42SBarry Smith pc[11] = m12 = p4*x9 + p8*x10 + p12*x11 + p16*x12; 8983287d42SBarry Smith 9083287d42SBarry Smith pc[12] = m13 = p1*x13 + p5*x14 + p9*x15 + p13*x16; 9183287d42SBarry Smith pc[13] = m14 = p2*x13 + p6*x14 + p10*x15 + p14*x16; 9283287d42SBarry Smith pc[14] = m15 = p3*x13 + p7*x14 + p11*x15 + p15*x16; 9383287d42SBarry Smith pc[15] = m16 = p4*x13 + p8*x14 + p12*x15 + p16*x16; 9483287d42SBarry Smith 9583287d42SBarry Smith nz = bi[row+1] - diag_offset[row] - 1; 9683287d42SBarry Smith pv += 16; 9783287d42SBarry Smith for (j=0; j<nz; j++) { 9883287d42SBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 9983287d42SBarry Smith x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 10083287d42SBarry Smith x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; 10183287d42SBarry Smith x14 = pv[13]; x15 = pv[14]; x16 = pv[15]; 10283287d42SBarry Smith x = rtmp + 16*pj[j]; 10383287d42SBarry Smith x[0] -= m1*x1 + m5*x2 + m9*x3 + m13*x4; 10483287d42SBarry Smith x[1] -= m2*x1 + m6*x2 + m10*x3 + m14*x4; 10583287d42SBarry Smith x[2] -= m3*x1 + m7*x2 + m11*x3 + m15*x4; 10683287d42SBarry Smith x[3] -= m4*x1 + m8*x2 + m12*x3 + m16*x4; 10783287d42SBarry Smith 10883287d42SBarry Smith x[4] -= m1*x5 + m5*x6 + m9*x7 + m13*x8; 10983287d42SBarry Smith x[5] -= m2*x5 + m6*x6 + m10*x7 + m14*x8; 11083287d42SBarry Smith x[6] -= m3*x5 + m7*x6 + m11*x7 + m15*x8; 11183287d42SBarry Smith x[7] -= m4*x5 + m8*x6 + m12*x7 + m16*x8; 11283287d42SBarry Smith 11383287d42SBarry Smith x[8] -= m1*x9 + m5*x10 + m9*x11 + m13*x12; 11483287d42SBarry Smith x[9] -= m2*x9 + m6*x10 + m10*x11 + m14*x12; 11583287d42SBarry Smith x[10] -= m3*x9 + m7*x10 + m11*x11 + m15*x12; 11683287d42SBarry Smith x[11] -= m4*x9 + m8*x10 + m12*x11 + m16*x12; 11783287d42SBarry Smith 11883287d42SBarry Smith x[12] -= m1*x13 + m5*x14 + m9*x15 + m13*x16; 11983287d42SBarry Smith x[13] -= m2*x13 + m6*x14 + m10*x15 + m14*x16; 12083287d42SBarry Smith x[14] -= m3*x13 + m7*x14 + m11*x15 + m15*x16; 12183287d42SBarry Smith x[15] -= m4*x13 + m8*x14 + m12*x15 + m16*x16; 12283287d42SBarry Smith 12383287d42SBarry Smith pv += 16; 12483287d42SBarry Smith } 125dc0b31edSSatish Balay ierr = PetscLogFlops(128.0*nz+112.0);CHKERRQ(ierr); 12683287d42SBarry Smith } 12783287d42SBarry Smith row = *ajtmp++; 12883287d42SBarry Smith } 12983287d42SBarry Smith /* finished row so stick it into b->a */ 13083287d42SBarry Smith pv = ba + 16*bi[i]; 13183287d42SBarry Smith pj = bj + bi[i]; 13283287d42SBarry Smith nz = bi[i+1] - bi[i]; 13383287d42SBarry Smith for (j=0; j<nz; j++) { 13483287d42SBarry Smith x = rtmp+16*pj[j]; 13583287d42SBarry Smith pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 13683287d42SBarry Smith pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8]; 13783287d42SBarry Smith pv[9] = x[9]; pv[10] = x[10]; pv[11] = x[11]; pv[12] = x[12]; 13883287d42SBarry Smith pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15]; 13983287d42SBarry Smith pv += 16; 14083287d42SBarry Smith } 14183287d42SBarry Smith /* invert diagonal block */ 14283287d42SBarry Smith w = ba + 16*diag_offset[i]; 143bcd9e38bSBarry Smith if (pivotinblocks) { 144*96b95a6bSBarry Smith ierr = PetscKernel_A_gets_inverse_A_4(w,shift);CHKERRQ(ierr); 145bcd9e38bSBarry Smith } else { 146*96b95a6bSBarry Smith ierr = PetscKernel_A_gets_inverse_A_4_nopivot(w);CHKERRQ(ierr); 147bcd9e38bSBarry Smith } 14883287d42SBarry Smith } 14983287d42SBarry Smith 15083287d42SBarry Smith ierr = PetscFree(rtmp);CHKERRQ(ierr); 15183287d42SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 15283287d42SBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 15306e38f1dSHong Zhang C->ops->solve = MatSolve_SeqBAIJ_4_inplace; 15406e38f1dSHong Zhang C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_4_inplace; 15583287d42SBarry Smith C->assembled = PETSC_TRUE; 156766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*4*4*4*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 15783287d42SBarry Smith PetscFunctionReturn(0); 15883287d42SBarry Smith } 159e6580ceeSShri Abhyankar 1604dd39f65SShri Abhyankar /* MatLUFactorNumeric_SeqBAIJ_4 - 1614dd39f65SShri Abhyankar copied from MatLUFactorNumeric_SeqBAIJ_N_inplace() and manually re-implemented 162*96b95a6bSBarry Smith PetscKernel_A_gets_A_times_B() 163*96b95a6bSBarry Smith PetscKernel_A_gets_A_minus_B_times_C() 164*96b95a6bSBarry Smith PetscKernel_A_gets_inverse_A() 165209027a4SShri Abhyankar */ 166c0c7eb62SShri Abhyankar 167209027a4SShri Abhyankar #undef __FUNCT__ 1684dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_4" 1694dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4(Mat B,Mat A,const MatFactorInfo *info) 170209027a4SShri Abhyankar { 171209027a4SShri Abhyankar Mat C=B; 172209027a4SShri Abhyankar Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; 173209027a4SShri Abhyankar IS isrow = b->row,isicol = b->icol; 174209027a4SShri Abhyankar PetscErrorCode ierr; 1755a586d82SBarry Smith const PetscInt *r,*ic; 176bbd65245SShri Abhyankar PetscInt i,j,k,nz,nzL,row; 177bbd65245SShri Abhyankar const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 178bbd65245SShri Abhyankar const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2; 179209027a4SShri Abhyankar MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 180bbd65245SShri Abhyankar PetscInt flg; 1816ba06ab7SHong Zhang PetscReal shift; 182209027a4SShri Abhyankar 183209027a4SShri Abhyankar PetscFunctionBegin; 184209027a4SShri Abhyankar ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 185209027a4SShri Abhyankar ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 186209027a4SShri Abhyankar 1876ba06ab7SHong Zhang if (info->shifttype == MAT_SHIFT_NONE){ 1886ba06ab7SHong Zhang shift = 0; 1896ba06ab7SHong Zhang } else { /* info->shifttype == MAT_SHIFT_INBLOCKS */ 1906ba06ab7SHong Zhang shift = info->shiftamount; 1916ba06ab7SHong Zhang } 1926ba06ab7SHong Zhang 193209027a4SShri Abhyankar /* generate work space needed by the factorization */ 194fca92195SBarry Smith ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr); 195209027a4SShri Abhyankar ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 196209027a4SShri Abhyankar 197209027a4SShri Abhyankar for (i=0; i<n; i++){ 198209027a4SShri Abhyankar /* zero rtmp */ 199209027a4SShri Abhyankar /* L part */ 200209027a4SShri Abhyankar nz = bi[i+1] - bi[i]; 201209027a4SShri Abhyankar bjtmp = bj + bi[i]; 202209027a4SShri Abhyankar for (j=0; j<nz; j++){ 203209027a4SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 204209027a4SShri Abhyankar } 205209027a4SShri Abhyankar 206209027a4SShri Abhyankar /* U part */ 20778bb4007SShri Abhyankar nz = bdiag[i]-bdiag[i+1]; 20878bb4007SShri Abhyankar bjtmp = bj + bdiag[i+1]+1; 20978bb4007SShri Abhyankar for (j=0; j<nz; j++){ 21078bb4007SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 21178bb4007SShri Abhyankar } 21278bb4007SShri Abhyankar 21378bb4007SShri Abhyankar /* load in initial (unfactored row) */ 21478bb4007SShri Abhyankar nz = ai[r[i]+1] - ai[r[i]]; 21578bb4007SShri Abhyankar ajtmp = aj + ai[r[i]]; 21678bb4007SShri Abhyankar v = aa + bs2*ai[r[i]]; 21778bb4007SShri Abhyankar for (j=0; j<nz; j++) { 21878bb4007SShri Abhyankar ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 21978bb4007SShri Abhyankar } 22078bb4007SShri Abhyankar 22178bb4007SShri Abhyankar /* elimination */ 22278bb4007SShri Abhyankar bjtmp = bj + bi[i]; 22378bb4007SShri Abhyankar nzL = bi[i+1] - bi[i]; 22478bb4007SShri Abhyankar for(k=0;k < nzL;k++) { 22578bb4007SShri Abhyankar row = bjtmp[k]; 22678bb4007SShri Abhyankar pc = rtmp + bs2*row; 22778bb4007SShri Abhyankar for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }} 22878bb4007SShri Abhyankar if (flg) { 22978bb4007SShri Abhyankar pv = b->a + bs2*bdiag[row]; 230*96b95a6bSBarry Smith /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 231*96b95a6bSBarry Smith ierr = PetscKernel_A_gets_A_times_B_4(pc,pv,mwork);CHKERRQ(ierr); 23278bb4007SShri Abhyankar 23378bb4007SShri Abhyankar pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */ 23478bb4007SShri Abhyankar pv = b->a + bs2*(bdiag[row+1]+1); 23578bb4007SShri Abhyankar nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */ 23678bb4007SShri Abhyankar for (j=0; j<nz; j++) { 237*96b95a6bSBarry Smith /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 23878bb4007SShri Abhyankar /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 23978bb4007SShri Abhyankar v = rtmp + bs2*pj[j]; 240*96b95a6bSBarry Smith ierr = PetscKernel_A_gets_A_minus_B_times_C_4(v,pc,pv);CHKERRQ(ierr); 24178bb4007SShri Abhyankar pv += bs2; 24278bb4007SShri Abhyankar } 24378bb4007SShri Abhyankar ierr = PetscLogFlops(128*nz+112);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 24478bb4007SShri Abhyankar } 24578bb4007SShri Abhyankar } 24678bb4007SShri Abhyankar 24778bb4007SShri Abhyankar /* finished row so stick it into b->a */ 24878bb4007SShri Abhyankar /* L part */ 24978bb4007SShri Abhyankar pv = b->a + bs2*bi[i] ; 25078bb4007SShri Abhyankar pj = b->j + bi[i] ; 25178bb4007SShri Abhyankar nz = bi[i+1] - bi[i]; 25278bb4007SShri Abhyankar for (j=0; j<nz; j++) { 25378bb4007SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 25478bb4007SShri Abhyankar } 25578bb4007SShri Abhyankar 25678bb4007SShri Abhyankar /* Mark diagonal and invert diagonal for simplier triangular solves */ 25778bb4007SShri Abhyankar pv = b->a + bs2*bdiag[i]; 25878bb4007SShri Abhyankar pj = b->j + bdiag[i]; 25978bb4007SShri Abhyankar ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 260*96b95a6bSBarry Smith ierr = PetscKernel_A_gets_inverse_A_4(pv,shift);CHKERRQ(ierr); 26178bb4007SShri Abhyankar 26278bb4007SShri Abhyankar /* U part */ 26378bb4007SShri Abhyankar pv = b->a + bs2*(bdiag[i+1]+1); 26478bb4007SShri Abhyankar pj = b->j + bdiag[i+1]+1; 26578bb4007SShri Abhyankar nz = bdiag[i] - bdiag[i+1] - 1; 26678bb4007SShri Abhyankar for (j=0; j<nz; j++){ 26778bb4007SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 26878bb4007SShri Abhyankar } 26978bb4007SShri Abhyankar } 27078bb4007SShri Abhyankar 271fca92195SBarry Smith ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); 27278bb4007SShri Abhyankar ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 27378bb4007SShri Abhyankar ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 2744dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_4; 2754dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_4; 27678bb4007SShri Abhyankar C->assembled = PETSC_TRUE; 277766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*4*4*4*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 27878bb4007SShri Abhyankar PetscFunctionReturn(0); 27978bb4007SShri Abhyankar } 28078bb4007SShri Abhyankar 28178bb4007SShri Abhyankar #undef __FUNCT__ 28206e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_inplace" 28306e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info) 284e6580ceeSShri Abhyankar { 285e6580ceeSShri Abhyankar /* 286e6580ceeSShri Abhyankar Default Version for when blocks are 4 by 4 Using natural ordering 287e6580ceeSShri Abhyankar */ 288e6580ceeSShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data; 289e6580ceeSShri Abhyankar PetscErrorCode ierr; 290e6580ceeSShri Abhyankar PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 291e6580ceeSShri Abhyankar PetscInt *ajtmpold,*ajtmp,nz,row; 292e6580ceeSShri Abhyankar PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; 293e6580ceeSShri Abhyankar MatScalar *pv,*v,*rtmp,*pc,*w,*x; 294e6580ceeSShri Abhyankar MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4; 295e6580ceeSShri Abhyankar MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,x16; 296e6580ceeSShri Abhyankar MatScalar p10,p11,p12,p13,p14,p15,p16,m10,m11,m12; 297e6580ceeSShri Abhyankar MatScalar m13,m14,m15,m16; 298e6580ceeSShri Abhyankar MatScalar *ba = b->a,*aa = a->a; 299ace3abfcSBarry Smith PetscBool pivotinblocks = b->pivotinblocks; 300182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 301e6580ceeSShri Abhyankar 302e6580ceeSShri Abhyankar PetscFunctionBegin; 303e6580ceeSShri Abhyankar ierr = PetscMalloc(16*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 304e6580ceeSShri Abhyankar 305e6580ceeSShri Abhyankar for (i=0; i<n; i++) { 306e6580ceeSShri Abhyankar nz = bi[i+1] - bi[i]; 307e6580ceeSShri Abhyankar ajtmp = bj + bi[i]; 308e6580ceeSShri Abhyankar for (j=0; j<nz; j++) { 309e6580ceeSShri Abhyankar x = rtmp+16*ajtmp[j]; 310e6580ceeSShri Abhyankar x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0; 311e6580ceeSShri Abhyankar x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = 0.0; 312e6580ceeSShri Abhyankar } 313e6580ceeSShri Abhyankar /* load in initial (unfactored row) */ 314e6580ceeSShri Abhyankar nz = ai[i+1] - ai[i]; 315e6580ceeSShri Abhyankar ajtmpold = aj + ai[i]; 316e6580ceeSShri Abhyankar v = aa + 16*ai[i]; 317e6580ceeSShri Abhyankar for (j=0; j<nz; j++) { 318e6580ceeSShri Abhyankar x = rtmp+16*ajtmpold[j]; 319e6580ceeSShri Abhyankar x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 320e6580ceeSShri Abhyankar x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8]; 321e6580ceeSShri Abhyankar x[9] = v[9]; x[10] = v[10]; x[11] = v[11]; x[12] = v[12]; x[13] = v[13]; 322e6580ceeSShri Abhyankar x[14] = v[14]; x[15] = v[15]; 323e6580ceeSShri Abhyankar v += 16; 324e6580ceeSShri Abhyankar } 325e6580ceeSShri Abhyankar row = *ajtmp++; 326e6580ceeSShri Abhyankar while (row < i) { 327e6580ceeSShri Abhyankar pc = rtmp + 16*row; 328e6580ceeSShri Abhyankar p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 329e6580ceeSShri Abhyankar p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8]; 330e6580ceeSShri Abhyankar p10 = pc[9]; p11 = pc[10]; p12 = pc[11]; p13 = pc[12]; p14 = pc[13]; 331e6580ceeSShri Abhyankar p15 = pc[14]; p16 = pc[15]; 332e6580ceeSShri Abhyankar if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || 333e6580ceeSShri Abhyankar p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 || 334e6580ceeSShri Abhyankar p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0 335e6580ceeSShri Abhyankar || p16 != 0.0) { 336e6580ceeSShri Abhyankar pv = ba + 16*diag_offset[row]; 337e6580ceeSShri Abhyankar pj = bj + diag_offset[row] + 1; 338e6580ceeSShri Abhyankar x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 339e6580ceeSShri Abhyankar x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 340e6580ceeSShri Abhyankar x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; x14 = pv[13]; 341e6580ceeSShri Abhyankar x15 = pv[14]; x16 = pv[15]; 342e6580ceeSShri Abhyankar pc[0] = m1 = p1*x1 + p5*x2 + p9*x3 + p13*x4; 343e6580ceeSShri Abhyankar pc[1] = m2 = p2*x1 + p6*x2 + p10*x3 + p14*x4; 344e6580ceeSShri Abhyankar pc[2] = m3 = p3*x1 + p7*x2 + p11*x3 + p15*x4; 345e6580ceeSShri Abhyankar pc[3] = m4 = p4*x1 + p8*x2 + p12*x3 + p16*x4; 346e6580ceeSShri Abhyankar 347e6580ceeSShri Abhyankar pc[4] = m5 = p1*x5 + p5*x6 + p9*x7 + p13*x8; 348e6580ceeSShri Abhyankar pc[5] = m6 = p2*x5 + p6*x6 + p10*x7 + p14*x8; 349e6580ceeSShri Abhyankar pc[6] = m7 = p3*x5 + p7*x6 + p11*x7 + p15*x8; 350e6580ceeSShri Abhyankar pc[7] = m8 = p4*x5 + p8*x6 + p12*x7 + p16*x8; 351e6580ceeSShri Abhyankar 352e6580ceeSShri Abhyankar pc[8] = m9 = p1*x9 + p5*x10 + p9*x11 + p13*x12; 353e6580ceeSShri Abhyankar pc[9] = m10 = p2*x9 + p6*x10 + p10*x11 + p14*x12; 354e6580ceeSShri Abhyankar pc[10] = m11 = p3*x9 + p7*x10 + p11*x11 + p15*x12; 355e6580ceeSShri Abhyankar pc[11] = m12 = p4*x9 + p8*x10 + p12*x11 + p16*x12; 356e6580ceeSShri Abhyankar 357e6580ceeSShri Abhyankar pc[12] = m13 = p1*x13 + p5*x14 + p9*x15 + p13*x16; 358e6580ceeSShri Abhyankar pc[13] = m14 = p2*x13 + p6*x14 + p10*x15 + p14*x16; 359e6580ceeSShri Abhyankar pc[14] = m15 = p3*x13 + p7*x14 + p11*x15 + p15*x16; 360e6580ceeSShri Abhyankar pc[15] = m16 = p4*x13 + p8*x14 + p12*x15 + p16*x16; 361e6580ceeSShri Abhyankar nz = bi[row+1] - diag_offset[row] - 1; 362e6580ceeSShri Abhyankar pv += 16; 363e6580ceeSShri Abhyankar for (j=0; j<nz; j++) { 364e6580ceeSShri Abhyankar x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 365e6580ceeSShri Abhyankar x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 366e6580ceeSShri Abhyankar x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; 367e6580ceeSShri Abhyankar x14 = pv[13]; x15 = pv[14]; x16 = pv[15]; 368e6580ceeSShri Abhyankar x = rtmp + 16*pj[j]; 369e6580ceeSShri Abhyankar x[0] -= m1*x1 + m5*x2 + m9*x3 + m13*x4; 370e6580ceeSShri Abhyankar x[1] -= m2*x1 + m6*x2 + m10*x3 + m14*x4; 371e6580ceeSShri Abhyankar x[2] -= m3*x1 + m7*x2 + m11*x3 + m15*x4; 372e6580ceeSShri Abhyankar x[3] -= m4*x1 + m8*x2 + m12*x3 + m16*x4; 373e6580ceeSShri Abhyankar 374e6580ceeSShri Abhyankar x[4] -= m1*x5 + m5*x6 + m9*x7 + m13*x8; 375e6580ceeSShri Abhyankar x[5] -= m2*x5 + m6*x6 + m10*x7 + m14*x8; 376e6580ceeSShri Abhyankar x[6] -= m3*x5 + m7*x6 + m11*x7 + m15*x8; 377e6580ceeSShri Abhyankar x[7] -= m4*x5 + m8*x6 + m12*x7 + m16*x8; 378e6580ceeSShri Abhyankar 379e6580ceeSShri Abhyankar x[8] -= m1*x9 + m5*x10 + m9*x11 + m13*x12; 380e6580ceeSShri Abhyankar x[9] -= m2*x9 + m6*x10 + m10*x11 + m14*x12; 381e6580ceeSShri Abhyankar x[10] -= m3*x9 + m7*x10 + m11*x11 + m15*x12; 382e6580ceeSShri Abhyankar x[11] -= m4*x9 + m8*x10 + m12*x11 + m16*x12; 383e6580ceeSShri Abhyankar 384e6580ceeSShri Abhyankar x[12] -= m1*x13 + m5*x14 + m9*x15 + m13*x16; 385e6580ceeSShri Abhyankar x[13] -= m2*x13 + m6*x14 + m10*x15 + m14*x16; 386e6580ceeSShri Abhyankar x[14] -= m3*x13 + m7*x14 + m11*x15 + m15*x16; 387e6580ceeSShri Abhyankar x[15] -= m4*x13 + m8*x14 + m12*x15 + m16*x16; 388e6580ceeSShri Abhyankar 389e6580ceeSShri Abhyankar pv += 16; 390e6580ceeSShri Abhyankar } 391e6580ceeSShri Abhyankar ierr = PetscLogFlops(128.0*nz+112.0);CHKERRQ(ierr); 392e6580ceeSShri Abhyankar } 393e6580ceeSShri Abhyankar row = *ajtmp++; 394e6580ceeSShri Abhyankar } 395e6580ceeSShri Abhyankar /* finished row so stick it into b->a */ 396e6580ceeSShri Abhyankar pv = ba + 16*bi[i]; 397e6580ceeSShri Abhyankar pj = bj + bi[i]; 398e6580ceeSShri Abhyankar nz = bi[i+1] - bi[i]; 399e6580ceeSShri Abhyankar for (j=0; j<nz; j++) { 400e6580ceeSShri Abhyankar x = rtmp+16*pj[j]; 401e6580ceeSShri Abhyankar pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 402e6580ceeSShri Abhyankar pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8]; 403e6580ceeSShri Abhyankar pv[9] = x[9]; pv[10] = x[10]; pv[11] = x[11]; pv[12] = x[12]; 404e6580ceeSShri Abhyankar pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15]; 405e6580ceeSShri Abhyankar pv += 16; 406e6580ceeSShri Abhyankar } 407e6580ceeSShri Abhyankar /* invert diagonal block */ 408e6580ceeSShri Abhyankar w = ba + 16*diag_offset[i]; 409e6580ceeSShri Abhyankar if (pivotinblocks) { 410*96b95a6bSBarry Smith ierr = PetscKernel_A_gets_inverse_A_4(w,shift);CHKERRQ(ierr); 411e6580ceeSShri Abhyankar } else { 412*96b95a6bSBarry Smith ierr = PetscKernel_A_gets_inverse_A_4_nopivot(w);CHKERRQ(ierr); 413e6580ceeSShri Abhyankar } 414e6580ceeSShri Abhyankar } 415e6580ceeSShri Abhyankar 416e6580ceeSShri Abhyankar ierr = PetscFree(rtmp);CHKERRQ(ierr); 41706e38f1dSHong Zhang C->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering_inplace; 41806e38f1dSHong Zhang C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_4_NaturalOrdering_inplace; 419e6580ceeSShri Abhyankar C->assembled = PETSC_TRUE; 420766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*4*4*4*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 421e6580ceeSShri Abhyankar PetscFunctionReturn(0); 422e6580ceeSShri Abhyankar } 423c0c7eb62SShri Abhyankar 424209027a4SShri Abhyankar /* 4254dd39f65SShri Abhyankar MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering - 4264dd39f65SShri Abhyankar copied from MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace() 427209027a4SShri Abhyankar */ 428209027a4SShri Abhyankar #undef __FUNCT__ 4294dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering" 4304dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info) 431209027a4SShri Abhyankar { 432209027a4SShri Abhyankar Mat C=B; 433209027a4SShri Abhyankar Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; 434209027a4SShri Abhyankar PetscErrorCode ierr; 435bbd65245SShri Abhyankar PetscInt i,j,k,nz,nzL,row; 436bbd65245SShri Abhyankar const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 437bbd65245SShri Abhyankar const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2; 438209027a4SShri Abhyankar MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 439bbd65245SShri Abhyankar PetscInt flg; 4406ba06ab7SHong Zhang PetscReal shift; 441e6580ceeSShri Abhyankar 442209027a4SShri Abhyankar PetscFunctionBegin; 443209027a4SShri Abhyankar /* generate work space needed by the factorization */ 444fca92195SBarry Smith ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr); 445209027a4SShri Abhyankar ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 446209027a4SShri Abhyankar 4476ba06ab7SHong Zhang if (info->shifttype == MAT_SHIFT_NONE){ 4486ba06ab7SHong Zhang shift = 0; 4496ba06ab7SHong Zhang } else { /* info->shifttype == MAT_SHIFT_INBLOCKS */ 4506ba06ab7SHong Zhang shift = info->shiftamount; 4516ba06ab7SHong Zhang } 4526ba06ab7SHong Zhang 453209027a4SShri Abhyankar for (i=0; i<n; i++){ 454209027a4SShri Abhyankar /* zero rtmp */ 455209027a4SShri Abhyankar /* L part */ 456209027a4SShri Abhyankar nz = bi[i+1] - bi[i]; 457209027a4SShri Abhyankar bjtmp = bj + bi[i]; 458209027a4SShri Abhyankar for (j=0; j<nz; j++){ 459209027a4SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 460209027a4SShri Abhyankar } 461209027a4SShri Abhyankar 462209027a4SShri Abhyankar /* U part */ 463b2b2dd24SShri Abhyankar nz = bdiag[i] - bdiag[i+1]; 464b2b2dd24SShri Abhyankar bjtmp = bj + bdiag[i+1]+1; 465b2b2dd24SShri Abhyankar for (j=0; j<nz; j++){ 466b2b2dd24SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 467b2b2dd24SShri Abhyankar } 468b2b2dd24SShri Abhyankar 469b2b2dd24SShri Abhyankar /* load in initial (unfactored row) */ 470b2b2dd24SShri Abhyankar nz = ai[i+1] - ai[i]; 471b2b2dd24SShri Abhyankar ajtmp = aj + ai[i]; 472b2b2dd24SShri Abhyankar v = aa + bs2*ai[i]; 473b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 474b2b2dd24SShri Abhyankar ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 475b2b2dd24SShri Abhyankar } 476b2b2dd24SShri Abhyankar 477b2b2dd24SShri Abhyankar /* elimination */ 478b2b2dd24SShri Abhyankar bjtmp = bj + bi[i]; 479b2b2dd24SShri Abhyankar nzL = bi[i+1] - bi[i]; 480b2b2dd24SShri Abhyankar for(k=0;k < nzL;k++) { 481b2b2dd24SShri Abhyankar row = bjtmp[k]; 482b2b2dd24SShri Abhyankar pc = rtmp + bs2*row; 483b2b2dd24SShri Abhyankar for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }} 484b2b2dd24SShri Abhyankar if (flg) { 485b2b2dd24SShri Abhyankar pv = b->a + bs2*bdiag[row]; 486*96b95a6bSBarry Smith /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 487*96b95a6bSBarry Smith ierr = PetscKernel_A_gets_A_times_B_4(pc,pv,mwork);CHKERRQ(ierr); 488b2b2dd24SShri Abhyankar 489b2b2dd24SShri Abhyankar pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */ 490b2b2dd24SShri Abhyankar pv = b->a + bs2*(bdiag[row+1]+1); 491b2b2dd24SShri Abhyankar nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */ 492b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 493*96b95a6bSBarry Smith /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 494b2b2dd24SShri Abhyankar /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 495b2b2dd24SShri Abhyankar v = rtmp + bs2*pj[j]; 496*96b95a6bSBarry Smith ierr = PetscKernel_A_gets_A_minus_B_times_C_4(v,pc,pv);CHKERRQ(ierr); 497b2b2dd24SShri Abhyankar pv += bs2; 498b2b2dd24SShri Abhyankar } 499b2b2dd24SShri Abhyankar ierr = PetscLogFlops(128*nz+112);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 500b2b2dd24SShri Abhyankar } 501b2b2dd24SShri Abhyankar } 502b2b2dd24SShri Abhyankar 503b2b2dd24SShri Abhyankar /* finished row so stick it into b->a */ 504b2b2dd24SShri Abhyankar /* L part */ 505b2b2dd24SShri Abhyankar pv = b->a + bs2*bi[i] ; 506b2b2dd24SShri Abhyankar pj = b->j + bi[i] ; 507b2b2dd24SShri Abhyankar nz = bi[i+1] - bi[i]; 508b2b2dd24SShri Abhyankar for (j=0; j<nz; j++) { 509b2b2dd24SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 510b2b2dd24SShri Abhyankar } 511b2b2dd24SShri Abhyankar 512b2b2dd24SShri Abhyankar /* Mark diagonal and invert diagonal for simplier triangular solves */ 513b2b2dd24SShri Abhyankar pv = b->a + bs2*bdiag[i]; 514b2b2dd24SShri Abhyankar pj = b->j + bdiag[i]; 515b2b2dd24SShri Abhyankar ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 516*96b95a6bSBarry Smith ierr = PetscKernel_A_gets_inverse_A_4(pv,shift);CHKERRQ(ierr); 517b2b2dd24SShri Abhyankar 518b2b2dd24SShri Abhyankar /* U part */ 519b2b2dd24SShri Abhyankar pv = b->a + bs2*(bdiag[i+1]+1); 520b2b2dd24SShri Abhyankar pj = b->j + bdiag[i+1]+1; 521b2b2dd24SShri Abhyankar nz = bdiag[i] - bdiag[i+1] - 1; 522b2b2dd24SShri Abhyankar for (j=0; j<nz; j++){ 523b2b2dd24SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 524b2b2dd24SShri Abhyankar } 525b2b2dd24SShri Abhyankar } 526fca92195SBarry Smith ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); 5274dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering; 5284dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_4_NaturalOrdering; 529b2b2dd24SShri Abhyankar C->assembled = PETSC_TRUE; 530766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*4*4*4*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 531b2b2dd24SShri Abhyankar PetscFunctionReturn(0); 532b2b2dd24SShri Abhyankar } 533b2b2dd24SShri Abhyankar 534e6580ceeSShri Abhyankar #if defined(PETSC_HAVE_SSE) 535e6580ceeSShri Abhyankar 536e6580ceeSShri Abhyankar #include PETSC_HAVE_SSE 537e6580ceeSShri Abhyankar 538e6580ceeSShri Abhyankar /* SSE Version for when blocks are 4 by 4 Using natural ordering */ 539e6580ceeSShri Abhyankar #undef __FUNCT__ 540e6580ceeSShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE" 541e6580ceeSShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE(Mat B,Mat A,const MatFactorInfo *info) 542e6580ceeSShri Abhyankar { 543e6580ceeSShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data; 544e6580ceeSShri Abhyankar PetscErrorCode ierr; 545e6580ceeSShri Abhyankar int i,j,n = a->mbs; 546e6580ceeSShri Abhyankar int *bj = b->j,*bjtmp,*pj; 547e6580ceeSShri Abhyankar int row; 548e6580ceeSShri Abhyankar int *ajtmpold,nz,*bi=b->i; 549e6580ceeSShri Abhyankar int *diag_offset = b->diag,*ai=a->i,*aj=a->j; 550e6580ceeSShri Abhyankar MatScalar *pv,*v,*rtmp,*pc,*w,*x; 551e6580ceeSShri Abhyankar MatScalar *ba = b->a,*aa = a->a; 552e6580ceeSShri Abhyankar int nonzero=0; 553e6580ceeSShri Abhyankar /* int nonzero=0,colscale = 16; */ 554ace3abfcSBarry Smith PetscBool pivotinblocks = b->pivotinblocks; 555182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 556e6580ceeSShri Abhyankar 557e6580ceeSShri Abhyankar PetscFunctionBegin; 558e6580ceeSShri Abhyankar SSE_SCOPE_BEGIN; 559e6580ceeSShri Abhyankar 560e32f2f54SBarry Smith if ((unsigned long)aa%16!=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Pointer aa is not 16 byte aligned. SSE will not work."); 561e32f2f54SBarry Smith if ((unsigned long)ba%16!=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Pointer ba is not 16 byte aligned. SSE will not work."); 562e6580ceeSShri Abhyankar ierr = PetscMalloc(16*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 563e32f2f54SBarry Smith if ((unsigned long)rtmp%16!=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Pointer rtmp is not 16 byte aligned. SSE will not work."); 564e6580ceeSShri Abhyankar /* if ((unsigned long)bj==(unsigned long)aj) { */ 565e6580ceeSShri Abhyankar /* colscale = 4; */ 566e6580ceeSShri Abhyankar /* } */ 567e6580ceeSShri Abhyankar for (i=0; i<n; i++) { 568e6580ceeSShri Abhyankar nz = bi[i+1] - bi[i]; 569e6580ceeSShri Abhyankar bjtmp = bj + bi[i]; 570e6580ceeSShri Abhyankar /* zero out the 4x4 block accumulators */ 571e6580ceeSShri Abhyankar /* zero out one register */ 572e6580ceeSShri Abhyankar XOR_PS(XMM7,XMM7); 573e6580ceeSShri Abhyankar for (j=0; j<nz; j++) { 574e6580ceeSShri Abhyankar x = rtmp+16*bjtmp[j]; 575e6580ceeSShri Abhyankar /* x = rtmp+4*bjtmp[j]; */ 576e6580ceeSShri Abhyankar SSE_INLINE_BEGIN_1(x) 577e6580ceeSShri Abhyankar /* Copy zero register to memory locations */ 578e6580ceeSShri Abhyankar /* Note: on future SSE architectures, STORE might be more efficient than STOREL/H */ 579e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_0,XMM7) 580e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_2,XMM7) 581e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_4,XMM7) 582e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_6,XMM7) 583e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_8,XMM7) 584e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_10,XMM7) 585e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_12,XMM7) 586e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_14,XMM7) 587e6580ceeSShri Abhyankar SSE_INLINE_END_1; 588e6580ceeSShri Abhyankar } 589e6580ceeSShri Abhyankar /* load in initial (unfactored row) */ 590e6580ceeSShri Abhyankar nz = ai[i+1] - ai[i]; 591e6580ceeSShri Abhyankar ajtmpold = aj + ai[i]; 592e6580ceeSShri Abhyankar v = aa + 16*ai[i]; 593e6580ceeSShri Abhyankar for (j=0; j<nz; j++) { 594e6580ceeSShri Abhyankar x = rtmp+16*ajtmpold[j]; 595e6580ceeSShri Abhyankar /* x = rtmp+colscale*ajtmpold[j]; */ 596e6580ceeSShri Abhyankar /* Copy v block into x block */ 597e6580ceeSShri Abhyankar SSE_INLINE_BEGIN_2(v,x) 598e6580ceeSShri Abhyankar /* Note: on future SSE architectures, STORE might be more efficient than STOREL/H */ 599e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_0,XMM0) 600e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_0,XMM0) 601e6580ceeSShri Abhyankar 602e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_2,XMM1) 603e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_2,XMM1) 604e6580ceeSShri Abhyankar 605e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_4,XMM2) 606e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_4,XMM2) 607e6580ceeSShri Abhyankar 608e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_6,XMM3) 609e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_6,XMM3) 610e6580ceeSShri Abhyankar 611e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_8,XMM4) 612e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_8,XMM4) 613e6580ceeSShri Abhyankar 614e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_10,XMM5) 615e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_10,XMM5) 616e6580ceeSShri Abhyankar 617e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_12,XMM6) 618e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_12,XMM6) 619e6580ceeSShri Abhyankar 620e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_14,XMM0) 621e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_14,XMM0) 622e6580ceeSShri Abhyankar SSE_INLINE_END_2; 623e6580ceeSShri Abhyankar 624e6580ceeSShri Abhyankar v += 16; 625e6580ceeSShri Abhyankar } 626e6580ceeSShri Abhyankar /* row = (*bjtmp++)/4; */ 627e6580ceeSShri Abhyankar row = *bjtmp++; 628e6580ceeSShri Abhyankar while (row < i) { 629e6580ceeSShri Abhyankar pc = rtmp + 16*row; 630e6580ceeSShri Abhyankar SSE_INLINE_BEGIN_1(pc) 631e6580ceeSShri Abhyankar /* Load block from lower triangle */ 632e6580ceeSShri Abhyankar /* Note: on future SSE architectures, STORE might be more efficient than STOREL/H */ 633e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_0,XMM0) 634e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_2,XMM0) 635e6580ceeSShri Abhyankar 636e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_4,XMM1) 637e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_6,XMM1) 638e6580ceeSShri Abhyankar 639e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_8,XMM2) 640e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_10,XMM2) 641e6580ceeSShri Abhyankar 642e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_12,XMM3) 643e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_14,XMM3) 644e6580ceeSShri Abhyankar 645e6580ceeSShri Abhyankar /* Compare block to zero block */ 646e6580ceeSShri Abhyankar 647e6580ceeSShri Abhyankar SSE_COPY_PS(XMM4,XMM7) 648e6580ceeSShri Abhyankar SSE_CMPNEQ_PS(XMM4,XMM0) 649e6580ceeSShri Abhyankar 650e6580ceeSShri Abhyankar SSE_COPY_PS(XMM5,XMM7) 651e6580ceeSShri Abhyankar SSE_CMPNEQ_PS(XMM5,XMM1) 652e6580ceeSShri Abhyankar 653e6580ceeSShri Abhyankar SSE_COPY_PS(XMM6,XMM7) 654e6580ceeSShri Abhyankar SSE_CMPNEQ_PS(XMM6,XMM2) 655e6580ceeSShri Abhyankar 656e6580ceeSShri Abhyankar SSE_CMPNEQ_PS(XMM7,XMM3) 657e6580ceeSShri Abhyankar 658e6580ceeSShri Abhyankar /* Reduce the comparisons to one SSE register */ 659e6580ceeSShri Abhyankar SSE_OR_PS(XMM6,XMM7) 660e6580ceeSShri Abhyankar SSE_OR_PS(XMM5,XMM4) 661e6580ceeSShri Abhyankar SSE_OR_PS(XMM5,XMM6) 662e6580ceeSShri Abhyankar SSE_INLINE_END_1; 663e6580ceeSShri Abhyankar 664e6580ceeSShri Abhyankar /* Reduce the one SSE register to an integer register for branching */ 665e6580ceeSShri Abhyankar /* Note: Since nonzero is an int, there is no INLINE block version of this call */ 666e6580ceeSShri Abhyankar MOVEMASK(nonzero,XMM5); 667e6580ceeSShri Abhyankar 668e6580ceeSShri Abhyankar /* If block is nonzero ... */ 669e6580ceeSShri Abhyankar if (nonzero) { 670e6580ceeSShri Abhyankar pv = ba + 16*diag_offset[row]; 671e6580ceeSShri Abhyankar PREFETCH_L1(&pv[16]); 672e6580ceeSShri Abhyankar pj = bj + diag_offset[row] + 1; 673e6580ceeSShri Abhyankar 674e6580ceeSShri Abhyankar /* Form Multiplier, one column at a time (Matrix-Matrix Product) */ 675e6580ceeSShri Abhyankar /* L_ij^(k+1) = L_ij^(k)*inv(L_jj^(k)) */ 676e6580ceeSShri Abhyankar /* but the diagonal was inverted already */ 677e6580ceeSShri Abhyankar /* and, L_ij^(k) is already loaded into registers XMM0-XMM3 columnwise */ 678e6580ceeSShri Abhyankar 679e6580ceeSShri Abhyankar SSE_INLINE_BEGIN_2(pv,pc) 680e6580ceeSShri Abhyankar /* Column 0, product is accumulated in XMM4 */ 681e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_0,XMM4) 682e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM4,XMM4,0x00) 683e6580ceeSShri Abhyankar SSE_MULT_PS(XMM4,XMM0) 684e6580ceeSShri Abhyankar 685e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_1,XMM5) 686e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 687e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM1) 688e6580ceeSShri Abhyankar SSE_ADD_PS(XMM4,XMM5) 689e6580ceeSShri Abhyankar 690e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_2,XMM6) 691e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 692e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM2) 693e6580ceeSShri Abhyankar SSE_ADD_PS(XMM4,XMM6) 694e6580ceeSShri Abhyankar 695e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_3,XMM7) 696e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 697e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM3) 698e6580ceeSShri Abhyankar SSE_ADD_PS(XMM4,XMM7) 699e6580ceeSShri Abhyankar 700e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_0,XMM4) 701e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_2,XMM4) 702e6580ceeSShri Abhyankar 703e6580ceeSShri Abhyankar /* Column 1, product is accumulated in XMM5 */ 704e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_4,XMM5) 705e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 706e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM0) 707e6580ceeSShri Abhyankar 708e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_5,XMM6) 709e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 710e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM1) 711e6580ceeSShri Abhyankar SSE_ADD_PS(XMM5,XMM6) 712e6580ceeSShri Abhyankar 713e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_6,XMM7) 714e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 715e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM2) 716e6580ceeSShri Abhyankar SSE_ADD_PS(XMM5,XMM7) 717e6580ceeSShri Abhyankar 718e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_7,XMM6) 719e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 720e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM3) 721e6580ceeSShri Abhyankar SSE_ADD_PS(XMM5,XMM6) 722e6580ceeSShri Abhyankar 723e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_4,XMM5) 724e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_6,XMM5) 725e6580ceeSShri Abhyankar 726e6580ceeSShri Abhyankar SSE_PREFETCH_L1(SSE_ARG_1,FLOAT_24) 727e6580ceeSShri Abhyankar 728e6580ceeSShri Abhyankar /* Column 2, product is accumulated in XMM6 */ 729e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_8,XMM6) 730e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 731e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM0) 732e6580ceeSShri Abhyankar 733e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_9,XMM7) 734e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 735e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM1) 736e6580ceeSShri Abhyankar SSE_ADD_PS(XMM6,XMM7) 737e6580ceeSShri Abhyankar 738e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_10,XMM7) 739e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 740e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM2) 741e6580ceeSShri Abhyankar SSE_ADD_PS(XMM6,XMM7) 742e6580ceeSShri Abhyankar 743e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_11,XMM7) 744e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 745e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM3) 746e6580ceeSShri Abhyankar SSE_ADD_PS(XMM6,XMM7) 747e6580ceeSShri Abhyankar 748e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_8,XMM6) 749e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_10,XMM6) 750e6580ceeSShri Abhyankar 751e6580ceeSShri Abhyankar /* Note: For the last column, we no longer need to preserve XMM0->XMM3 */ 752e6580ceeSShri Abhyankar /* Column 3, product is accumulated in XMM0 */ 753e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_12,XMM7) 754e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 755e6580ceeSShri Abhyankar SSE_MULT_PS(XMM0,XMM7) 756e6580ceeSShri Abhyankar 757e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_13,XMM7) 758e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 759e6580ceeSShri Abhyankar SSE_MULT_PS(XMM1,XMM7) 760e6580ceeSShri Abhyankar SSE_ADD_PS(XMM0,XMM1) 761e6580ceeSShri Abhyankar 762e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_14,XMM1) 763e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM1,XMM1,0x00) 764e6580ceeSShri Abhyankar SSE_MULT_PS(XMM1,XMM2) 765e6580ceeSShri Abhyankar SSE_ADD_PS(XMM0,XMM1) 766e6580ceeSShri Abhyankar 767e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_15,XMM7) 768e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 769e6580ceeSShri Abhyankar SSE_MULT_PS(XMM3,XMM7) 770e6580ceeSShri Abhyankar SSE_ADD_PS(XMM0,XMM3) 771e6580ceeSShri Abhyankar 772e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_12,XMM0) 773e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_14,XMM0) 774e6580ceeSShri Abhyankar 775e6580ceeSShri Abhyankar /* Simplify Bookkeeping -- Completely Unnecessary Instructions */ 776e6580ceeSShri Abhyankar /* This is code to be maintained and read by humans afterall. */ 777e6580ceeSShri Abhyankar /* Copy Multiplier Col 3 into XMM3 */ 778e6580ceeSShri Abhyankar SSE_COPY_PS(XMM3,XMM0) 779e6580ceeSShri Abhyankar /* Copy Multiplier Col 2 into XMM2 */ 780e6580ceeSShri Abhyankar SSE_COPY_PS(XMM2,XMM6) 781e6580ceeSShri Abhyankar /* Copy Multiplier Col 1 into XMM1 */ 782e6580ceeSShri Abhyankar SSE_COPY_PS(XMM1,XMM5) 783e6580ceeSShri Abhyankar /* Copy Multiplier Col 0 into XMM0 */ 784e6580ceeSShri Abhyankar SSE_COPY_PS(XMM0,XMM4) 785e6580ceeSShri Abhyankar SSE_INLINE_END_2; 786e6580ceeSShri Abhyankar 787e6580ceeSShri Abhyankar /* Update the row: */ 788e6580ceeSShri Abhyankar nz = bi[row+1] - diag_offset[row] - 1; 789e6580ceeSShri Abhyankar pv += 16; 790e6580ceeSShri Abhyankar for (j=0; j<nz; j++) { 791e6580ceeSShri Abhyankar PREFETCH_L1(&pv[16]); 792e6580ceeSShri Abhyankar x = rtmp + 16*pj[j]; 793e6580ceeSShri Abhyankar /* x = rtmp + 4*pj[j]; */ 794e6580ceeSShri Abhyankar 795e6580ceeSShri Abhyankar /* X:=X-M*PV, One column at a time */ 796e6580ceeSShri Abhyankar /* Note: M is already loaded columnwise into registers XMM0-XMM3 */ 797e6580ceeSShri Abhyankar SSE_INLINE_BEGIN_2(x,pv) 798e6580ceeSShri Abhyankar /* Load First Column of X*/ 799e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_0,XMM4) 800e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_2,XMM4) 801e6580ceeSShri Abhyankar 802e6580ceeSShri Abhyankar /* Matrix-Vector Product: */ 803e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_0,XMM5) 804e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 805e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM0) 806e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM5) 807e6580ceeSShri Abhyankar 808e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_1,XMM6) 809e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 810e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM1) 811e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM6) 812e6580ceeSShri Abhyankar 813e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_2,XMM7) 814e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 815e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM2) 816e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM7) 817e6580ceeSShri Abhyankar 818e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_3,XMM5) 819e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 820e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM3) 821e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM5) 822e6580ceeSShri Abhyankar 823e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_0,XMM4) 824e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_2,XMM4) 825e6580ceeSShri Abhyankar 826e6580ceeSShri Abhyankar /* Second Column */ 827e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_4,XMM5) 828e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_6,XMM5) 829e6580ceeSShri Abhyankar 830e6580ceeSShri Abhyankar /* Matrix-Vector Product: */ 831e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_4,XMM6) 832e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 833e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM0) 834e6580ceeSShri Abhyankar SSE_SUB_PS(XMM5,XMM6) 835e6580ceeSShri Abhyankar 836e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_5,XMM7) 837e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 838e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM1) 839e6580ceeSShri Abhyankar SSE_SUB_PS(XMM5,XMM7) 840e6580ceeSShri Abhyankar 841e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_6,XMM4) 842e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM4,XMM4,0x00) 843e6580ceeSShri Abhyankar SSE_MULT_PS(XMM4,XMM2) 844e6580ceeSShri Abhyankar SSE_SUB_PS(XMM5,XMM4) 845e6580ceeSShri Abhyankar 846e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_7,XMM6) 847e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 848e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM3) 849e6580ceeSShri Abhyankar SSE_SUB_PS(XMM5,XMM6) 850e6580ceeSShri Abhyankar 851e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_4,XMM5) 852e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_6,XMM5) 853e6580ceeSShri Abhyankar 854e6580ceeSShri Abhyankar SSE_PREFETCH_L1(SSE_ARG_2,FLOAT_24) 855e6580ceeSShri Abhyankar 856e6580ceeSShri Abhyankar /* Third Column */ 857e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_8,XMM6) 858e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_10,XMM6) 859e6580ceeSShri Abhyankar 860e6580ceeSShri Abhyankar /* Matrix-Vector Product: */ 861e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_8,XMM7) 862e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 863e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM0) 864e6580ceeSShri Abhyankar SSE_SUB_PS(XMM6,XMM7) 865e6580ceeSShri Abhyankar 866e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_9,XMM4) 867e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM4,XMM4,0x00) 868e6580ceeSShri Abhyankar SSE_MULT_PS(XMM4,XMM1) 869e6580ceeSShri Abhyankar SSE_SUB_PS(XMM6,XMM4) 870e6580ceeSShri Abhyankar 871e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_10,XMM5) 872e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 873e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM2) 874e6580ceeSShri Abhyankar SSE_SUB_PS(XMM6,XMM5) 875e6580ceeSShri Abhyankar 876e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_11,XMM7) 877e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 878e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM3) 879e6580ceeSShri Abhyankar SSE_SUB_PS(XMM6,XMM7) 880e6580ceeSShri Abhyankar 881e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_8,XMM6) 882e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_10,XMM6) 883e6580ceeSShri Abhyankar 884e6580ceeSShri Abhyankar /* Fourth Column */ 885e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_12,XMM4) 886e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_14,XMM4) 887e6580ceeSShri Abhyankar 888e6580ceeSShri Abhyankar /* Matrix-Vector Product: */ 889e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_12,XMM5) 890e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 891e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM0) 892e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM5) 893e6580ceeSShri Abhyankar 894e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_13,XMM6) 895e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 896e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM1) 897e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM6) 898e6580ceeSShri Abhyankar 899e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_14,XMM7) 900e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 901e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM2) 902e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM7) 903e6580ceeSShri Abhyankar 904e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_15,XMM5) 905e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 906e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM3) 907e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM5) 908e6580ceeSShri Abhyankar 909e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_12,XMM4) 910e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_14,XMM4) 911e6580ceeSShri Abhyankar SSE_INLINE_END_2; 912e6580ceeSShri Abhyankar pv += 16; 913e6580ceeSShri Abhyankar } 914e6580ceeSShri Abhyankar ierr = PetscLogFlops(128.0*nz+112.0);CHKERRQ(ierr); 915e6580ceeSShri Abhyankar } 916e6580ceeSShri Abhyankar row = *bjtmp++; 917e6580ceeSShri Abhyankar /* row = (*bjtmp++)/4; */ 918e6580ceeSShri Abhyankar } 919e6580ceeSShri Abhyankar /* finished row so stick it into b->a */ 920e6580ceeSShri Abhyankar pv = ba + 16*bi[i]; 921e6580ceeSShri Abhyankar pj = bj + bi[i]; 922e6580ceeSShri Abhyankar nz = bi[i+1] - bi[i]; 923e6580ceeSShri Abhyankar 924e6580ceeSShri Abhyankar /* Copy x block back into pv block */ 925e6580ceeSShri Abhyankar for (j=0; j<nz; j++) { 926e6580ceeSShri Abhyankar x = rtmp+16*pj[j]; 927e6580ceeSShri Abhyankar /* x = rtmp+4*pj[j]; */ 928e6580ceeSShri Abhyankar 929e6580ceeSShri Abhyankar SSE_INLINE_BEGIN_2(x,pv) 930e6580ceeSShri Abhyankar /* Note: on future SSE architectures, STORE might be more efficient than STOREL/H */ 931e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_0,XMM1) 932e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_0,XMM1) 933e6580ceeSShri Abhyankar 934e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_2,XMM2) 935e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_2,XMM2) 936e6580ceeSShri Abhyankar 937e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_4,XMM3) 938e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_4,XMM3) 939e6580ceeSShri Abhyankar 940e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_6,XMM4) 941e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_6,XMM4) 942e6580ceeSShri Abhyankar 943e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_8,XMM5) 944e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_8,XMM5) 945e6580ceeSShri Abhyankar 946e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_10,XMM6) 947e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_10,XMM6) 948e6580ceeSShri Abhyankar 949e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_12,XMM7) 950e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_12,XMM7) 951e6580ceeSShri Abhyankar 952e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_14,XMM0) 953e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_14,XMM0) 954e6580ceeSShri Abhyankar SSE_INLINE_END_2; 955e6580ceeSShri Abhyankar pv += 16; 956e6580ceeSShri Abhyankar } 957e6580ceeSShri Abhyankar /* invert diagonal block */ 958e6580ceeSShri Abhyankar w = ba + 16*diag_offset[i]; 959e6580ceeSShri Abhyankar if (pivotinblocks) { 960*96b95a6bSBarry Smith ierr = PetscKernel_A_gets_inverse_A_4(w,shift);CHKERRQ(ierr); 961e6580ceeSShri Abhyankar } else { 962*96b95a6bSBarry Smith ierr = PetscKernel_A_gets_inverse_A_4_nopivot(w);CHKERRQ(ierr); 963e6580ceeSShri Abhyankar } 964*96b95a6bSBarry Smith /* ierr = PetscKernel_A_gets_inverse_A_4_SSE(w);CHKERRQ(ierr); */ 965e6580ceeSShri Abhyankar /* Note: Using Kramer's rule, flop count below might be infairly high or low? */ 966e6580ceeSShri Abhyankar } 967e6580ceeSShri Abhyankar 968e6580ceeSShri Abhyankar ierr = PetscFree(rtmp);CHKERRQ(ierr); 969e6580ceeSShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering_SSE; 970e6580ceeSShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_4_NaturalOrdering_SSE; 971e6580ceeSShri Abhyankar C->assembled = PETSC_TRUE; 972766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*bs*bs2*b->mbs);CHKERRQ(ierr); 973e6580ceeSShri Abhyankar /* Flop Count from inverting diagonal blocks */ 974e6580ceeSShri Abhyankar SSE_SCOPE_END; 975e6580ceeSShri Abhyankar PetscFunctionReturn(0); 976e6580ceeSShri Abhyankar } 977e6580ceeSShri Abhyankar 978e6580ceeSShri Abhyankar #undef __FUNCT__ 979e6580ceeSShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj_Inplace" 980e6580ceeSShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj_Inplace(Mat C) 981e6580ceeSShri Abhyankar { 982e6580ceeSShri Abhyankar Mat A=C; 983e6580ceeSShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data; 984e6580ceeSShri Abhyankar PetscErrorCode ierr; 985e6580ceeSShri Abhyankar int i,j,n = a->mbs; 986e6580ceeSShri Abhyankar unsigned short *bj = (unsigned short *)(b->j),*bjtmp,*pj; 987e6580ceeSShri Abhyankar unsigned short *aj = (unsigned short *)(a->j),*ajtmp; 988e6580ceeSShri Abhyankar unsigned int row; 989e6580ceeSShri Abhyankar int nz,*bi=b->i; 990e6580ceeSShri Abhyankar int *diag_offset = b->diag,*ai=a->i; 991e6580ceeSShri Abhyankar MatScalar *pv,*v,*rtmp,*pc,*w,*x; 992e6580ceeSShri Abhyankar MatScalar *ba = b->a,*aa = a->a; 993e6580ceeSShri Abhyankar int nonzero=0; 994e6580ceeSShri Abhyankar /* int nonzero=0,colscale = 16; */ 995ace3abfcSBarry Smith PetscBool pivotinblocks = b->pivotinblocks; 996182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 997e6580ceeSShri Abhyankar 998e6580ceeSShri Abhyankar PetscFunctionBegin; 999e6580ceeSShri Abhyankar SSE_SCOPE_BEGIN; 1000e6580ceeSShri Abhyankar 1001e32f2f54SBarry Smith if ((unsigned long)aa%16!=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Pointer aa is not 16 byte aligned. SSE will not work."); 1002e32f2f54SBarry Smith if ((unsigned long)ba%16!=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Pointer ba is not 16 byte aligned. SSE will not work."); 1003e6580ceeSShri Abhyankar ierr = PetscMalloc(16*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 1004e32f2f54SBarry Smith if ((unsigned long)rtmp%16!=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Pointer rtmp is not 16 byte aligned. SSE will not work."); 1005e6580ceeSShri Abhyankar /* if ((unsigned long)bj==(unsigned long)aj) { */ 1006e6580ceeSShri Abhyankar /* colscale = 4; */ 1007e6580ceeSShri Abhyankar /* } */ 1008e6580ceeSShri Abhyankar 1009e6580ceeSShri Abhyankar for (i=0; i<n; i++) { 1010e6580ceeSShri Abhyankar nz = bi[i+1] - bi[i]; 1011e6580ceeSShri Abhyankar bjtmp = bj + bi[i]; 1012e6580ceeSShri Abhyankar /* zero out the 4x4 block accumulators */ 1013e6580ceeSShri Abhyankar /* zero out one register */ 1014e6580ceeSShri Abhyankar XOR_PS(XMM7,XMM7); 1015e6580ceeSShri Abhyankar for (j=0; j<nz; j++) { 1016e6580ceeSShri Abhyankar x = rtmp+16*((unsigned int)bjtmp[j]); 1017e6580ceeSShri Abhyankar /* x = rtmp+4*bjtmp[j]; */ 1018e6580ceeSShri Abhyankar SSE_INLINE_BEGIN_1(x) 1019e6580ceeSShri Abhyankar /* Copy zero register to memory locations */ 1020e6580ceeSShri Abhyankar /* Note: on future SSE architectures, STORE might be more efficient than STOREL/H */ 1021e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_0,XMM7) 1022e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_2,XMM7) 1023e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_4,XMM7) 1024e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_6,XMM7) 1025e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_8,XMM7) 1026e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_10,XMM7) 1027e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_12,XMM7) 1028e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_14,XMM7) 1029e6580ceeSShri Abhyankar SSE_INLINE_END_1; 1030e6580ceeSShri Abhyankar } 1031e6580ceeSShri Abhyankar /* load in initial (unfactored row) */ 1032e6580ceeSShri Abhyankar nz = ai[i+1] - ai[i]; 1033e6580ceeSShri Abhyankar ajtmp = aj + ai[i]; 1034e6580ceeSShri Abhyankar v = aa + 16*ai[i]; 1035e6580ceeSShri Abhyankar for (j=0; j<nz; j++) { 1036e6580ceeSShri Abhyankar x = rtmp+16*((unsigned int)ajtmp[j]); 1037e6580ceeSShri Abhyankar /* x = rtmp+colscale*ajtmp[j]; */ 1038e6580ceeSShri Abhyankar /* Copy v block into x block */ 1039e6580ceeSShri Abhyankar SSE_INLINE_BEGIN_2(v,x) 1040e6580ceeSShri Abhyankar /* Note: on future SSE architectures, STORE might be more efficient than STOREL/H */ 1041e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_0,XMM0) 1042e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_0,XMM0) 1043e6580ceeSShri Abhyankar 1044e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_2,XMM1) 1045e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_2,XMM1) 1046e6580ceeSShri Abhyankar 1047e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_4,XMM2) 1048e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_4,XMM2) 1049e6580ceeSShri Abhyankar 1050e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_6,XMM3) 1051e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_6,XMM3) 1052e6580ceeSShri Abhyankar 1053e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_8,XMM4) 1054e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_8,XMM4) 1055e6580ceeSShri Abhyankar 1056e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_10,XMM5) 1057e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_10,XMM5) 1058e6580ceeSShri Abhyankar 1059e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_12,XMM6) 1060e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_12,XMM6) 1061e6580ceeSShri Abhyankar 1062e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_14,XMM0) 1063e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_14,XMM0) 1064e6580ceeSShri Abhyankar SSE_INLINE_END_2; 1065e6580ceeSShri Abhyankar 1066e6580ceeSShri Abhyankar v += 16; 1067e6580ceeSShri Abhyankar } 1068e6580ceeSShri Abhyankar /* row = (*bjtmp++)/4; */ 1069e6580ceeSShri Abhyankar row = (unsigned int)(*bjtmp++); 1070e6580ceeSShri Abhyankar while (row < i) { 1071e6580ceeSShri Abhyankar pc = rtmp + 16*row; 1072e6580ceeSShri Abhyankar SSE_INLINE_BEGIN_1(pc) 1073e6580ceeSShri Abhyankar /* Load block from lower triangle */ 1074e6580ceeSShri Abhyankar /* Note: on future SSE architectures, STORE might be more efficient than STOREL/H */ 1075e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_0,XMM0) 1076e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_2,XMM0) 1077e6580ceeSShri Abhyankar 1078e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_4,XMM1) 1079e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_6,XMM1) 1080e6580ceeSShri Abhyankar 1081e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_8,XMM2) 1082e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_10,XMM2) 1083e6580ceeSShri Abhyankar 1084e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_12,XMM3) 1085e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_14,XMM3) 1086e6580ceeSShri Abhyankar 1087e6580ceeSShri Abhyankar /* Compare block to zero block */ 1088e6580ceeSShri Abhyankar 1089e6580ceeSShri Abhyankar SSE_COPY_PS(XMM4,XMM7) 1090e6580ceeSShri Abhyankar SSE_CMPNEQ_PS(XMM4,XMM0) 1091e6580ceeSShri Abhyankar 1092e6580ceeSShri Abhyankar SSE_COPY_PS(XMM5,XMM7) 1093e6580ceeSShri Abhyankar SSE_CMPNEQ_PS(XMM5,XMM1) 1094e6580ceeSShri Abhyankar 1095e6580ceeSShri Abhyankar SSE_COPY_PS(XMM6,XMM7) 1096e6580ceeSShri Abhyankar SSE_CMPNEQ_PS(XMM6,XMM2) 1097e6580ceeSShri Abhyankar 1098e6580ceeSShri Abhyankar SSE_CMPNEQ_PS(XMM7,XMM3) 1099e6580ceeSShri Abhyankar 1100e6580ceeSShri Abhyankar /* Reduce the comparisons to one SSE register */ 1101e6580ceeSShri Abhyankar SSE_OR_PS(XMM6,XMM7) 1102e6580ceeSShri Abhyankar SSE_OR_PS(XMM5,XMM4) 1103e6580ceeSShri Abhyankar SSE_OR_PS(XMM5,XMM6) 1104e6580ceeSShri Abhyankar SSE_INLINE_END_1; 1105e6580ceeSShri Abhyankar 1106e6580ceeSShri Abhyankar /* Reduce the one SSE register to an integer register for branching */ 1107e6580ceeSShri Abhyankar /* Note: Since nonzero is an int, there is no INLINE block version of this call */ 1108e6580ceeSShri Abhyankar MOVEMASK(nonzero,XMM5); 1109e6580ceeSShri Abhyankar 1110e6580ceeSShri Abhyankar /* If block is nonzero ... */ 1111e6580ceeSShri Abhyankar if (nonzero) { 1112e6580ceeSShri Abhyankar pv = ba + 16*diag_offset[row]; 1113e6580ceeSShri Abhyankar PREFETCH_L1(&pv[16]); 1114e6580ceeSShri Abhyankar pj = bj + diag_offset[row] + 1; 1115e6580ceeSShri Abhyankar 1116e6580ceeSShri Abhyankar /* Form Multiplier, one column at a time (Matrix-Matrix Product) */ 1117e6580ceeSShri Abhyankar /* L_ij^(k+1) = L_ij^(k)*inv(L_jj^(k)) */ 1118e6580ceeSShri Abhyankar /* but the diagonal was inverted already */ 1119e6580ceeSShri Abhyankar /* and, L_ij^(k) is already loaded into registers XMM0-XMM3 columnwise */ 1120e6580ceeSShri Abhyankar 1121e6580ceeSShri Abhyankar SSE_INLINE_BEGIN_2(pv,pc) 1122e6580ceeSShri Abhyankar /* Column 0, product is accumulated in XMM4 */ 1123e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_0,XMM4) 1124e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM4,XMM4,0x00) 1125e6580ceeSShri Abhyankar SSE_MULT_PS(XMM4,XMM0) 1126e6580ceeSShri Abhyankar 1127e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_1,XMM5) 1128e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 1129e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM1) 1130e6580ceeSShri Abhyankar SSE_ADD_PS(XMM4,XMM5) 1131e6580ceeSShri Abhyankar 1132e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_2,XMM6) 1133e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 1134e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM2) 1135e6580ceeSShri Abhyankar SSE_ADD_PS(XMM4,XMM6) 1136e6580ceeSShri Abhyankar 1137e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_3,XMM7) 1138e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1139e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM3) 1140e6580ceeSShri Abhyankar SSE_ADD_PS(XMM4,XMM7) 1141e6580ceeSShri Abhyankar 1142e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_0,XMM4) 1143e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_2,XMM4) 1144e6580ceeSShri Abhyankar 1145e6580ceeSShri Abhyankar /* Column 1, product is accumulated in XMM5 */ 1146e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_4,XMM5) 1147e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 1148e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM0) 1149e6580ceeSShri Abhyankar 1150e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_5,XMM6) 1151e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 1152e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM1) 1153e6580ceeSShri Abhyankar SSE_ADD_PS(XMM5,XMM6) 1154e6580ceeSShri Abhyankar 1155e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_6,XMM7) 1156e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1157e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM2) 1158e6580ceeSShri Abhyankar SSE_ADD_PS(XMM5,XMM7) 1159e6580ceeSShri Abhyankar 1160e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_7,XMM6) 1161e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 1162e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM3) 1163e6580ceeSShri Abhyankar SSE_ADD_PS(XMM5,XMM6) 1164e6580ceeSShri Abhyankar 1165e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_4,XMM5) 1166e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_6,XMM5) 1167e6580ceeSShri Abhyankar 1168e6580ceeSShri Abhyankar SSE_PREFETCH_L1(SSE_ARG_1,FLOAT_24) 1169e6580ceeSShri Abhyankar 1170e6580ceeSShri Abhyankar /* Column 2, product is accumulated in XMM6 */ 1171e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_8,XMM6) 1172e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 1173e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM0) 1174e6580ceeSShri Abhyankar 1175e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_9,XMM7) 1176e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1177e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM1) 1178e6580ceeSShri Abhyankar SSE_ADD_PS(XMM6,XMM7) 1179e6580ceeSShri Abhyankar 1180e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_10,XMM7) 1181e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1182e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM2) 1183e6580ceeSShri Abhyankar SSE_ADD_PS(XMM6,XMM7) 1184e6580ceeSShri Abhyankar 1185e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_11,XMM7) 1186e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1187e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM3) 1188e6580ceeSShri Abhyankar SSE_ADD_PS(XMM6,XMM7) 1189e6580ceeSShri Abhyankar 1190e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_8,XMM6) 1191e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_10,XMM6) 1192e6580ceeSShri Abhyankar 1193e6580ceeSShri Abhyankar /* Note: For the last column, we no longer need to preserve XMM0->XMM3 */ 1194e6580ceeSShri Abhyankar /* Column 3, product is accumulated in XMM0 */ 1195e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_12,XMM7) 1196e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1197e6580ceeSShri Abhyankar SSE_MULT_PS(XMM0,XMM7) 1198e6580ceeSShri Abhyankar 1199e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_13,XMM7) 1200e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1201e6580ceeSShri Abhyankar SSE_MULT_PS(XMM1,XMM7) 1202e6580ceeSShri Abhyankar SSE_ADD_PS(XMM0,XMM1) 1203e6580ceeSShri Abhyankar 1204e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_14,XMM1) 1205e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM1,XMM1,0x00) 1206e6580ceeSShri Abhyankar SSE_MULT_PS(XMM1,XMM2) 1207e6580ceeSShri Abhyankar SSE_ADD_PS(XMM0,XMM1) 1208e6580ceeSShri Abhyankar 1209e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_15,XMM7) 1210e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1211e6580ceeSShri Abhyankar SSE_MULT_PS(XMM3,XMM7) 1212e6580ceeSShri Abhyankar SSE_ADD_PS(XMM0,XMM3) 1213e6580ceeSShri Abhyankar 1214e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_12,XMM0) 1215e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_14,XMM0) 1216e6580ceeSShri Abhyankar 1217e6580ceeSShri Abhyankar /* Simplify Bookkeeping -- Completely Unnecessary Instructions */ 1218e6580ceeSShri Abhyankar /* This is code to be maintained and read by humans afterall. */ 1219e6580ceeSShri Abhyankar /* Copy Multiplier Col 3 into XMM3 */ 1220e6580ceeSShri Abhyankar SSE_COPY_PS(XMM3,XMM0) 1221e6580ceeSShri Abhyankar /* Copy Multiplier Col 2 into XMM2 */ 1222e6580ceeSShri Abhyankar SSE_COPY_PS(XMM2,XMM6) 1223e6580ceeSShri Abhyankar /* Copy Multiplier Col 1 into XMM1 */ 1224e6580ceeSShri Abhyankar SSE_COPY_PS(XMM1,XMM5) 1225e6580ceeSShri Abhyankar /* Copy Multiplier Col 0 into XMM0 */ 1226e6580ceeSShri Abhyankar SSE_COPY_PS(XMM0,XMM4) 1227e6580ceeSShri Abhyankar SSE_INLINE_END_2; 1228e6580ceeSShri Abhyankar 1229e6580ceeSShri Abhyankar /* Update the row: */ 1230e6580ceeSShri Abhyankar nz = bi[row+1] - diag_offset[row] - 1; 1231e6580ceeSShri Abhyankar pv += 16; 1232e6580ceeSShri Abhyankar for (j=0; j<nz; j++) { 1233e6580ceeSShri Abhyankar PREFETCH_L1(&pv[16]); 1234e6580ceeSShri Abhyankar x = rtmp + 16*((unsigned int)pj[j]); 1235e6580ceeSShri Abhyankar /* x = rtmp + 4*pj[j]; */ 1236e6580ceeSShri Abhyankar 1237e6580ceeSShri Abhyankar /* X:=X-M*PV, One column at a time */ 1238e6580ceeSShri Abhyankar /* Note: M is already loaded columnwise into registers XMM0-XMM3 */ 1239e6580ceeSShri Abhyankar SSE_INLINE_BEGIN_2(x,pv) 1240e6580ceeSShri Abhyankar /* Load First Column of X*/ 1241e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_0,XMM4) 1242e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_2,XMM4) 1243e6580ceeSShri Abhyankar 1244e6580ceeSShri Abhyankar /* Matrix-Vector Product: */ 1245e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_0,XMM5) 1246e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 1247e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM0) 1248e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM5) 1249e6580ceeSShri Abhyankar 1250e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_1,XMM6) 1251e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 1252e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM1) 1253e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM6) 1254e6580ceeSShri Abhyankar 1255e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_2,XMM7) 1256e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1257e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM2) 1258e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM7) 1259e6580ceeSShri Abhyankar 1260e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_3,XMM5) 1261e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 1262e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM3) 1263e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM5) 1264e6580ceeSShri Abhyankar 1265e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_0,XMM4) 1266e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_2,XMM4) 1267e6580ceeSShri Abhyankar 1268e6580ceeSShri Abhyankar /* Second Column */ 1269e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_4,XMM5) 1270e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_6,XMM5) 1271e6580ceeSShri Abhyankar 1272e6580ceeSShri Abhyankar /* Matrix-Vector Product: */ 1273e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_4,XMM6) 1274e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 1275e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM0) 1276e6580ceeSShri Abhyankar SSE_SUB_PS(XMM5,XMM6) 1277e6580ceeSShri Abhyankar 1278e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_5,XMM7) 1279e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1280e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM1) 1281e6580ceeSShri Abhyankar SSE_SUB_PS(XMM5,XMM7) 1282e6580ceeSShri Abhyankar 1283e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_6,XMM4) 1284e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM4,XMM4,0x00) 1285e6580ceeSShri Abhyankar SSE_MULT_PS(XMM4,XMM2) 1286e6580ceeSShri Abhyankar SSE_SUB_PS(XMM5,XMM4) 1287e6580ceeSShri Abhyankar 1288e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_7,XMM6) 1289e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 1290e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM3) 1291e6580ceeSShri Abhyankar SSE_SUB_PS(XMM5,XMM6) 1292e6580ceeSShri Abhyankar 1293e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_4,XMM5) 1294e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_6,XMM5) 1295e6580ceeSShri Abhyankar 1296e6580ceeSShri Abhyankar SSE_PREFETCH_L1(SSE_ARG_2,FLOAT_24) 1297e6580ceeSShri Abhyankar 1298e6580ceeSShri Abhyankar /* Third Column */ 1299e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_8,XMM6) 1300e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_10,XMM6) 1301e6580ceeSShri Abhyankar 1302e6580ceeSShri Abhyankar /* Matrix-Vector Product: */ 1303e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_8,XMM7) 1304e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1305e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM0) 1306e6580ceeSShri Abhyankar SSE_SUB_PS(XMM6,XMM7) 1307e6580ceeSShri Abhyankar 1308e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_9,XMM4) 1309e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM4,XMM4,0x00) 1310e6580ceeSShri Abhyankar SSE_MULT_PS(XMM4,XMM1) 1311e6580ceeSShri Abhyankar SSE_SUB_PS(XMM6,XMM4) 1312e6580ceeSShri Abhyankar 1313e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_10,XMM5) 1314e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 1315e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM2) 1316e6580ceeSShri Abhyankar SSE_SUB_PS(XMM6,XMM5) 1317e6580ceeSShri Abhyankar 1318e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_11,XMM7) 1319e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1320e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM3) 1321e6580ceeSShri Abhyankar SSE_SUB_PS(XMM6,XMM7) 1322e6580ceeSShri Abhyankar 1323e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_8,XMM6) 1324e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_10,XMM6) 1325e6580ceeSShri Abhyankar 1326e6580ceeSShri Abhyankar /* Fourth Column */ 1327e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_12,XMM4) 1328e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_14,XMM4) 1329e6580ceeSShri Abhyankar 1330e6580ceeSShri Abhyankar /* Matrix-Vector Product: */ 1331e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_12,XMM5) 1332e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 1333e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM0) 1334e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM5) 1335e6580ceeSShri Abhyankar 1336e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_13,XMM6) 1337e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 1338e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM1) 1339e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM6) 1340e6580ceeSShri Abhyankar 1341e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_14,XMM7) 1342e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1343e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM2) 1344e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM7) 1345e6580ceeSShri Abhyankar 1346e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_15,XMM5) 1347e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 1348e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM3) 1349e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM5) 1350e6580ceeSShri Abhyankar 1351e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_12,XMM4) 1352e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_14,XMM4) 1353e6580ceeSShri Abhyankar SSE_INLINE_END_2; 1354e6580ceeSShri Abhyankar pv += 16; 1355e6580ceeSShri Abhyankar } 1356e6580ceeSShri Abhyankar ierr = PetscLogFlops(128.0*nz+112.0);CHKERRQ(ierr); 1357e6580ceeSShri Abhyankar } 1358e6580ceeSShri Abhyankar row = (unsigned int)(*bjtmp++); 1359e6580ceeSShri Abhyankar /* row = (*bjtmp++)/4; */ 1360e6580ceeSShri Abhyankar /* bjtmp++; */ 1361e6580ceeSShri Abhyankar } 1362e6580ceeSShri Abhyankar /* finished row so stick it into b->a */ 1363e6580ceeSShri Abhyankar pv = ba + 16*bi[i]; 1364e6580ceeSShri Abhyankar pj = bj + bi[i]; 1365e6580ceeSShri Abhyankar nz = bi[i+1] - bi[i]; 1366e6580ceeSShri Abhyankar 1367e6580ceeSShri Abhyankar /* Copy x block back into pv block */ 1368e6580ceeSShri Abhyankar for (j=0; j<nz; j++) { 1369e6580ceeSShri Abhyankar x = rtmp+16*((unsigned int)pj[j]); 1370e6580ceeSShri Abhyankar /* x = rtmp+4*pj[j]; */ 1371e6580ceeSShri Abhyankar 1372e6580ceeSShri Abhyankar SSE_INLINE_BEGIN_2(x,pv) 1373e6580ceeSShri Abhyankar /* Note: on future SSE architectures, STORE might be more efficient than STOREL/H */ 1374e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_0,XMM1) 1375e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_0,XMM1) 1376e6580ceeSShri Abhyankar 1377e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_2,XMM2) 1378e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_2,XMM2) 1379e6580ceeSShri Abhyankar 1380e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_4,XMM3) 1381e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_4,XMM3) 1382e6580ceeSShri Abhyankar 1383e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_6,XMM4) 1384e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_6,XMM4) 1385e6580ceeSShri Abhyankar 1386e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_8,XMM5) 1387e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_8,XMM5) 1388e6580ceeSShri Abhyankar 1389e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_10,XMM6) 1390e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_10,XMM6) 1391e6580ceeSShri Abhyankar 1392e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_12,XMM7) 1393e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_12,XMM7) 1394e6580ceeSShri Abhyankar 1395e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_14,XMM0) 1396e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_14,XMM0) 1397e6580ceeSShri Abhyankar SSE_INLINE_END_2; 1398e6580ceeSShri Abhyankar pv += 16; 1399e6580ceeSShri Abhyankar } 1400e6580ceeSShri Abhyankar /* invert diagonal block */ 1401e6580ceeSShri Abhyankar w = ba + 16*diag_offset[i]; 1402e6580ceeSShri Abhyankar if (pivotinblocks) { 1403*96b95a6bSBarry Smith ierr = PetscKernel_A_gets_inverse_A_4(w,shift);CHKERRQ(ierr); 1404e6580ceeSShri Abhyankar } else { 1405*96b95a6bSBarry Smith ierr = PetscKernel_A_gets_inverse_A_4_nopivot(w);CHKERRQ(ierr); 1406e6580ceeSShri Abhyankar } 1407*96b95a6bSBarry Smith /* ierr = PetscKernel_A_gets_inverse_A_4_SSE(w);CHKERRQ(ierr); */ 1408e6580ceeSShri Abhyankar /* Note: Using Kramer's rule, flop count below might be infairly high or low? */ 1409e6580ceeSShri Abhyankar } 1410e6580ceeSShri Abhyankar 1411e6580ceeSShri Abhyankar ierr = PetscFree(rtmp);CHKERRQ(ierr); 1412e6580ceeSShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering_SSE; 1413e6580ceeSShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_4_NaturalOrdering_SSE; 1414e6580ceeSShri Abhyankar C->assembled = PETSC_TRUE; 1415766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*bs*bs2*b->mbs);CHKERRQ(ierr); 1416e6580ceeSShri Abhyankar /* Flop Count from inverting diagonal blocks */ 1417e6580ceeSShri Abhyankar SSE_SCOPE_END; 1418e6580ceeSShri Abhyankar PetscFunctionReturn(0); 1419e6580ceeSShri Abhyankar } 1420e6580ceeSShri Abhyankar 1421e6580ceeSShri Abhyankar #undef __FUNCT__ 1422e6580ceeSShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj" 1423e6580ceeSShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat C,Mat A,const MatFactorInfo *info) 1424e6580ceeSShri Abhyankar { 1425e6580ceeSShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data; 1426e6580ceeSShri Abhyankar PetscErrorCode ierr; 1427e6580ceeSShri Abhyankar int i,j,n = a->mbs; 1428e6580ceeSShri Abhyankar unsigned short *bj = (unsigned short *)(b->j),*bjtmp,*pj; 1429e6580ceeSShri Abhyankar unsigned int row; 1430e6580ceeSShri Abhyankar int *ajtmpold,nz,*bi=b->i; 1431e6580ceeSShri Abhyankar int *diag_offset = b->diag,*ai=a->i,*aj=a->j; 1432e6580ceeSShri Abhyankar MatScalar *pv,*v,*rtmp,*pc,*w,*x; 1433e6580ceeSShri Abhyankar MatScalar *ba = b->a,*aa = a->a; 1434e6580ceeSShri Abhyankar int nonzero=0; 1435e6580ceeSShri Abhyankar /* int nonzero=0,colscale = 16; */ 1436ace3abfcSBarry Smith PetscBool pivotinblocks = b->pivotinblocks; 1437182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 1438e6580ceeSShri Abhyankar 1439e6580ceeSShri Abhyankar PetscFunctionBegin; 1440e6580ceeSShri Abhyankar SSE_SCOPE_BEGIN; 1441e6580ceeSShri Abhyankar 1442e32f2f54SBarry Smith if ((unsigned long)aa%16!=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Pointer aa is not 16 byte aligned. SSE will not work."); 1443e32f2f54SBarry Smith if ((unsigned long)ba%16!=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Pointer ba is not 16 byte aligned. SSE will not work."); 1444e6580ceeSShri Abhyankar ierr = PetscMalloc(16*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 1445e32f2f54SBarry Smith if ((unsigned long)rtmp%16!=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Pointer rtmp is not 16 byte aligned. SSE will not work."); 1446e6580ceeSShri Abhyankar /* if ((unsigned long)bj==(unsigned long)aj) { */ 1447e6580ceeSShri Abhyankar /* colscale = 4; */ 1448e6580ceeSShri Abhyankar /* } */ 1449e6580ceeSShri Abhyankar if ((unsigned long)bj==(unsigned long)aj) { 1450e6580ceeSShri Abhyankar return(MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj_Inplace(C)); 1451e6580ceeSShri Abhyankar } 1452e6580ceeSShri Abhyankar 1453e6580ceeSShri Abhyankar for (i=0; i<n; i++) { 1454e6580ceeSShri Abhyankar nz = bi[i+1] - bi[i]; 1455e6580ceeSShri Abhyankar bjtmp = bj + bi[i]; 1456e6580ceeSShri Abhyankar /* zero out the 4x4 block accumulators */ 1457e6580ceeSShri Abhyankar /* zero out one register */ 1458e6580ceeSShri Abhyankar XOR_PS(XMM7,XMM7); 1459e6580ceeSShri Abhyankar for (j=0; j<nz; j++) { 1460e6580ceeSShri Abhyankar x = rtmp+16*((unsigned int)bjtmp[j]); 1461e6580ceeSShri Abhyankar /* x = rtmp+4*bjtmp[j]; */ 1462e6580ceeSShri Abhyankar SSE_INLINE_BEGIN_1(x) 1463e6580ceeSShri Abhyankar /* Copy zero register to memory locations */ 1464e6580ceeSShri Abhyankar /* Note: on future SSE architectures, STORE might be more efficient than STOREL/H */ 1465e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_0,XMM7) 1466e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_2,XMM7) 1467e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_4,XMM7) 1468e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_6,XMM7) 1469e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_8,XMM7) 1470e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_10,XMM7) 1471e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_12,XMM7) 1472e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_14,XMM7) 1473e6580ceeSShri Abhyankar SSE_INLINE_END_1; 1474e6580ceeSShri Abhyankar } 1475e6580ceeSShri Abhyankar /* load in initial (unfactored row) */ 1476e6580ceeSShri Abhyankar nz = ai[i+1] - ai[i]; 1477e6580ceeSShri Abhyankar ajtmpold = aj + ai[i]; 1478e6580ceeSShri Abhyankar v = aa + 16*ai[i]; 1479e6580ceeSShri Abhyankar for (j=0; j<nz; j++) { 1480e6580ceeSShri Abhyankar x = rtmp+16*ajtmpold[j]; 1481e6580ceeSShri Abhyankar /* x = rtmp+colscale*ajtmpold[j]; */ 1482e6580ceeSShri Abhyankar /* Copy v block into x block */ 1483e6580ceeSShri Abhyankar SSE_INLINE_BEGIN_2(v,x) 1484e6580ceeSShri Abhyankar /* Note: on future SSE architectures, STORE might be more efficient than STOREL/H */ 1485e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_0,XMM0) 1486e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_0,XMM0) 1487e6580ceeSShri Abhyankar 1488e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_2,XMM1) 1489e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_2,XMM1) 1490e6580ceeSShri Abhyankar 1491e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_4,XMM2) 1492e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_4,XMM2) 1493e6580ceeSShri Abhyankar 1494e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_6,XMM3) 1495e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_6,XMM3) 1496e6580ceeSShri Abhyankar 1497e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_8,XMM4) 1498e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_8,XMM4) 1499e6580ceeSShri Abhyankar 1500e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_10,XMM5) 1501e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_10,XMM5) 1502e6580ceeSShri Abhyankar 1503e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_12,XMM6) 1504e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_12,XMM6) 1505e6580ceeSShri Abhyankar 1506e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_14,XMM0) 1507e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_14,XMM0) 1508e6580ceeSShri Abhyankar SSE_INLINE_END_2; 1509e6580ceeSShri Abhyankar 1510e6580ceeSShri Abhyankar v += 16; 1511e6580ceeSShri Abhyankar } 1512e6580ceeSShri Abhyankar /* row = (*bjtmp++)/4; */ 1513e6580ceeSShri Abhyankar row = (unsigned int)(*bjtmp++); 1514e6580ceeSShri Abhyankar while (row < i) { 1515e6580ceeSShri Abhyankar pc = rtmp + 16*row; 1516e6580ceeSShri Abhyankar SSE_INLINE_BEGIN_1(pc) 1517e6580ceeSShri Abhyankar /* Load block from lower triangle */ 1518e6580ceeSShri Abhyankar /* Note: on future SSE architectures, STORE might be more efficient than STOREL/H */ 1519e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_0,XMM0) 1520e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_2,XMM0) 1521e6580ceeSShri Abhyankar 1522e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_4,XMM1) 1523e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_6,XMM1) 1524e6580ceeSShri Abhyankar 1525e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_8,XMM2) 1526e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_10,XMM2) 1527e6580ceeSShri Abhyankar 1528e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_12,XMM3) 1529e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_14,XMM3) 1530e6580ceeSShri Abhyankar 1531e6580ceeSShri Abhyankar /* Compare block to zero block */ 1532e6580ceeSShri Abhyankar 1533e6580ceeSShri Abhyankar SSE_COPY_PS(XMM4,XMM7) 1534e6580ceeSShri Abhyankar SSE_CMPNEQ_PS(XMM4,XMM0) 1535e6580ceeSShri Abhyankar 1536e6580ceeSShri Abhyankar SSE_COPY_PS(XMM5,XMM7) 1537e6580ceeSShri Abhyankar SSE_CMPNEQ_PS(XMM5,XMM1) 1538e6580ceeSShri Abhyankar 1539e6580ceeSShri Abhyankar SSE_COPY_PS(XMM6,XMM7) 1540e6580ceeSShri Abhyankar SSE_CMPNEQ_PS(XMM6,XMM2) 1541e6580ceeSShri Abhyankar 1542e6580ceeSShri Abhyankar SSE_CMPNEQ_PS(XMM7,XMM3) 1543e6580ceeSShri Abhyankar 1544e6580ceeSShri Abhyankar /* Reduce the comparisons to one SSE register */ 1545e6580ceeSShri Abhyankar SSE_OR_PS(XMM6,XMM7) 1546e6580ceeSShri Abhyankar SSE_OR_PS(XMM5,XMM4) 1547e6580ceeSShri Abhyankar SSE_OR_PS(XMM5,XMM6) 1548e6580ceeSShri Abhyankar SSE_INLINE_END_1; 1549e6580ceeSShri Abhyankar 1550e6580ceeSShri Abhyankar /* Reduce the one SSE register to an integer register for branching */ 1551e6580ceeSShri Abhyankar /* Note: Since nonzero is an int, there is no INLINE block version of this call */ 1552e6580ceeSShri Abhyankar MOVEMASK(nonzero,XMM5); 1553e6580ceeSShri Abhyankar 1554e6580ceeSShri Abhyankar /* If block is nonzero ... */ 1555e6580ceeSShri Abhyankar if (nonzero) { 1556e6580ceeSShri Abhyankar pv = ba + 16*diag_offset[row]; 1557e6580ceeSShri Abhyankar PREFETCH_L1(&pv[16]); 1558e6580ceeSShri Abhyankar pj = bj + diag_offset[row] + 1; 1559e6580ceeSShri Abhyankar 1560e6580ceeSShri Abhyankar /* Form Multiplier, one column at a time (Matrix-Matrix Product) */ 1561e6580ceeSShri Abhyankar /* L_ij^(k+1) = L_ij^(k)*inv(L_jj^(k)) */ 1562e6580ceeSShri Abhyankar /* but the diagonal was inverted already */ 1563e6580ceeSShri Abhyankar /* and, L_ij^(k) is already loaded into registers XMM0-XMM3 columnwise */ 1564e6580ceeSShri Abhyankar 1565e6580ceeSShri Abhyankar SSE_INLINE_BEGIN_2(pv,pc) 1566e6580ceeSShri Abhyankar /* Column 0, product is accumulated in XMM4 */ 1567e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_0,XMM4) 1568e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM4,XMM4,0x00) 1569e6580ceeSShri Abhyankar SSE_MULT_PS(XMM4,XMM0) 1570e6580ceeSShri Abhyankar 1571e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_1,XMM5) 1572e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 1573e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM1) 1574e6580ceeSShri Abhyankar SSE_ADD_PS(XMM4,XMM5) 1575e6580ceeSShri Abhyankar 1576e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_2,XMM6) 1577e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 1578e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM2) 1579e6580ceeSShri Abhyankar SSE_ADD_PS(XMM4,XMM6) 1580e6580ceeSShri Abhyankar 1581e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_3,XMM7) 1582e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1583e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM3) 1584e6580ceeSShri Abhyankar SSE_ADD_PS(XMM4,XMM7) 1585e6580ceeSShri Abhyankar 1586e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_0,XMM4) 1587e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_2,XMM4) 1588e6580ceeSShri Abhyankar 1589e6580ceeSShri Abhyankar /* Column 1, product is accumulated in XMM5 */ 1590e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_4,XMM5) 1591e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 1592e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM0) 1593e6580ceeSShri Abhyankar 1594e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_5,XMM6) 1595e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 1596e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM1) 1597e6580ceeSShri Abhyankar SSE_ADD_PS(XMM5,XMM6) 1598e6580ceeSShri Abhyankar 1599e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_6,XMM7) 1600e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1601e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM2) 1602e6580ceeSShri Abhyankar SSE_ADD_PS(XMM5,XMM7) 1603e6580ceeSShri Abhyankar 1604e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_7,XMM6) 1605e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 1606e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM3) 1607e6580ceeSShri Abhyankar SSE_ADD_PS(XMM5,XMM6) 1608e6580ceeSShri Abhyankar 1609e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_4,XMM5) 1610e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_6,XMM5) 1611e6580ceeSShri Abhyankar 1612e6580ceeSShri Abhyankar SSE_PREFETCH_L1(SSE_ARG_1,FLOAT_24) 1613e6580ceeSShri Abhyankar 1614e6580ceeSShri Abhyankar /* Column 2, product is accumulated in XMM6 */ 1615e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_8,XMM6) 1616e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 1617e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM0) 1618e6580ceeSShri Abhyankar 1619e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_9,XMM7) 1620e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1621e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM1) 1622e6580ceeSShri Abhyankar SSE_ADD_PS(XMM6,XMM7) 1623e6580ceeSShri Abhyankar 1624e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_10,XMM7) 1625e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1626e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM2) 1627e6580ceeSShri Abhyankar SSE_ADD_PS(XMM6,XMM7) 1628e6580ceeSShri Abhyankar 1629e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_11,XMM7) 1630e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1631e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM3) 1632e6580ceeSShri Abhyankar SSE_ADD_PS(XMM6,XMM7) 1633e6580ceeSShri Abhyankar 1634e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_8,XMM6) 1635e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_10,XMM6) 1636e6580ceeSShri Abhyankar 1637e6580ceeSShri Abhyankar /* Note: For the last column, we no longer need to preserve XMM0->XMM3 */ 1638e6580ceeSShri Abhyankar /* Column 3, product is accumulated in XMM0 */ 1639e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_12,XMM7) 1640e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1641e6580ceeSShri Abhyankar SSE_MULT_PS(XMM0,XMM7) 1642e6580ceeSShri Abhyankar 1643e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_13,XMM7) 1644e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1645e6580ceeSShri Abhyankar SSE_MULT_PS(XMM1,XMM7) 1646e6580ceeSShri Abhyankar SSE_ADD_PS(XMM0,XMM1) 1647e6580ceeSShri Abhyankar 1648e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_14,XMM1) 1649e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM1,XMM1,0x00) 1650e6580ceeSShri Abhyankar SSE_MULT_PS(XMM1,XMM2) 1651e6580ceeSShri Abhyankar SSE_ADD_PS(XMM0,XMM1) 1652e6580ceeSShri Abhyankar 1653e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_1,FLOAT_15,XMM7) 1654e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1655e6580ceeSShri Abhyankar SSE_MULT_PS(XMM3,XMM7) 1656e6580ceeSShri Abhyankar SSE_ADD_PS(XMM0,XMM3) 1657e6580ceeSShri Abhyankar 1658e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_12,XMM0) 1659e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_14,XMM0) 1660e6580ceeSShri Abhyankar 1661e6580ceeSShri Abhyankar /* Simplify Bookkeeping -- Completely Unnecessary Instructions */ 1662e6580ceeSShri Abhyankar /* This is code to be maintained and read by humans afterall. */ 1663e6580ceeSShri Abhyankar /* Copy Multiplier Col 3 into XMM3 */ 1664e6580ceeSShri Abhyankar SSE_COPY_PS(XMM3,XMM0) 1665e6580ceeSShri Abhyankar /* Copy Multiplier Col 2 into XMM2 */ 1666e6580ceeSShri Abhyankar SSE_COPY_PS(XMM2,XMM6) 1667e6580ceeSShri Abhyankar /* Copy Multiplier Col 1 into XMM1 */ 1668e6580ceeSShri Abhyankar SSE_COPY_PS(XMM1,XMM5) 1669e6580ceeSShri Abhyankar /* Copy Multiplier Col 0 into XMM0 */ 1670e6580ceeSShri Abhyankar SSE_COPY_PS(XMM0,XMM4) 1671e6580ceeSShri Abhyankar SSE_INLINE_END_2; 1672e6580ceeSShri Abhyankar 1673e6580ceeSShri Abhyankar /* Update the row: */ 1674e6580ceeSShri Abhyankar nz = bi[row+1] - diag_offset[row] - 1; 1675e6580ceeSShri Abhyankar pv += 16; 1676e6580ceeSShri Abhyankar for (j=0; j<nz; j++) { 1677e6580ceeSShri Abhyankar PREFETCH_L1(&pv[16]); 1678e6580ceeSShri Abhyankar x = rtmp + 16*((unsigned int)pj[j]); 1679e6580ceeSShri Abhyankar /* x = rtmp + 4*pj[j]; */ 1680e6580ceeSShri Abhyankar 1681e6580ceeSShri Abhyankar /* X:=X-M*PV, One column at a time */ 1682e6580ceeSShri Abhyankar /* Note: M is already loaded columnwise into registers XMM0-XMM3 */ 1683e6580ceeSShri Abhyankar SSE_INLINE_BEGIN_2(x,pv) 1684e6580ceeSShri Abhyankar /* Load First Column of X*/ 1685e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_0,XMM4) 1686e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_2,XMM4) 1687e6580ceeSShri Abhyankar 1688e6580ceeSShri Abhyankar /* Matrix-Vector Product: */ 1689e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_0,XMM5) 1690e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 1691e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM0) 1692e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM5) 1693e6580ceeSShri Abhyankar 1694e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_1,XMM6) 1695e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 1696e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM1) 1697e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM6) 1698e6580ceeSShri Abhyankar 1699e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_2,XMM7) 1700e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1701e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM2) 1702e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM7) 1703e6580ceeSShri Abhyankar 1704e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_3,XMM5) 1705e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 1706e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM3) 1707e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM5) 1708e6580ceeSShri Abhyankar 1709e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_0,XMM4) 1710e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_2,XMM4) 1711e6580ceeSShri Abhyankar 1712e6580ceeSShri Abhyankar /* Second Column */ 1713e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_4,XMM5) 1714e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_6,XMM5) 1715e6580ceeSShri Abhyankar 1716e6580ceeSShri Abhyankar /* Matrix-Vector Product: */ 1717e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_4,XMM6) 1718e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 1719e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM0) 1720e6580ceeSShri Abhyankar SSE_SUB_PS(XMM5,XMM6) 1721e6580ceeSShri Abhyankar 1722e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_5,XMM7) 1723e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1724e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM1) 1725e6580ceeSShri Abhyankar SSE_SUB_PS(XMM5,XMM7) 1726e6580ceeSShri Abhyankar 1727e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_6,XMM4) 1728e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM4,XMM4,0x00) 1729e6580ceeSShri Abhyankar SSE_MULT_PS(XMM4,XMM2) 1730e6580ceeSShri Abhyankar SSE_SUB_PS(XMM5,XMM4) 1731e6580ceeSShri Abhyankar 1732e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_7,XMM6) 1733e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 1734e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM3) 1735e6580ceeSShri Abhyankar SSE_SUB_PS(XMM5,XMM6) 1736e6580ceeSShri Abhyankar 1737e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_4,XMM5) 1738e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_6,XMM5) 1739e6580ceeSShri Abhyankar 1740e6580ceeSShri Abhyankar SSE_PREFETCH_L1(SSE_ARG_2,FLOAT_24) 1741e6580ceeSShri Abhyankar 1742e6580ceeSShri Abhyankar /* Third Column */ 1743e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_8,XMM6) 1744e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_10,XMM6) 1745e6580ceeSShri Abhyankar 1746e6580ceeSShri Abhyankar /* Matrix-Vector Product: */ 1747e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_8,XMM7) 1748e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1749e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM0) 1750e6580ceeSShri Abhyankar SSE_SUB_PS(XMM6,XMM7) 1751e6580ceeSShri Abhyankar 1752e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_9,XMM4) 1753e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM4,XMM4,0x00) 1754e6580ceeSShri Abhyankar SSE_MULT_PS(XMM4,XMM1) 1755e6580ceeSShri Abhyankar SSE_SUB_PS(XMM6,XMM4) 1756e6580ceeSShri Abhyankar 1757e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_10,XMM5) 1758e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 1759e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM2) 1760e6580ceeSShri Abhyankar SSE_SUB_PS(XMM6,XMM5) 1761e6580ceeSShri Abhyankar 1762e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_11,XMM7) 1763e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1764e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM3) 1765e6580ceeSShri Abhyankar SSE_SUB_PS(XMM6,XMM7) 1766e6580ceeSShri Abhyankar 1767e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_8,XMM6) 1768e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_10,XMM6) 1769e6580ceeSShri Abhyankar 1770e6580ceeSShri Abhyankar /* Fourth Column */ 1771e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_12,XMM4) 1772e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_14,XMM4) 1773e6580ceeSShri Abhyankar 1774e6580ceeSShri Abhyankar /* Matrix-Vector Product: */ 1775e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_12,XMM5) 1776e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 1777e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM0) 1778e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM5) 1779e6580ceeSShri Abhyankar 1780e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_13,XMM6) 1781e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM6,XMM6,0x00) 1782e6580ceeSShri Abhyankar SSE_MULT_PS(XMM6,XMM1) 1783e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM6) 1784e6580ceeSShri Abhyankar 1785e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_14,XMM7) 1786e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM7,XMM7,0x00) 1787e6580ceeSShri Abhyankar SSE_MULT_PS(XMM7,XMM2) 1788e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM7) 1789e6580ceeSShri Abhyankar 1790e6580ceeSShri Abhyankar SSE_LOAD_SS(SSE_ARG_2,FLOAT_15,XMM5) 1791e6580ceeSShri Abhyankar SSE_SHUFFLE(XMM5,XMM5,0x00) 1792e6580ceeSShri Abhyankar SSE_MULT_PS(XMM5,XMM3) 1793e6580ceeSShri Abhyankar SSE_SUB_PS(XMM4,XMM5) 1794e6580ceeSShri Abhyankar 1795e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_1,FLOAT_12,XMM4) 1796e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_1,FLOAT_14,XMM4) 1797e6580ceeSShri Abhyankar SSE_INLINE_END_2; 1798e6580ceeSShri Abhyankar pv += 16; 1799e6580ceeSShri Abhyankar } 1800e6580ceeSShri Abhyankar ierr = PetscLogFlops(128.0*nz+112.0);CHKERRQ(ierr); 1801e6580ceeSShri Abhyankar } 1802e6580ceeSShri Abhyankar row = (unsigned int)(*bjtmp++); 1803e6580ceeSShri Abhyankar /* row = (*bjtmp++)/4; */ 1804e6580ceeSShri Abhyankar /* bjtmp++; */ 1805e6580ceeSShri Abhyankar } 1806e6580ceeSShri Abhyankar /* finished row so stick it into b->a */ 1807e6580ceeSShri Abhyankar pv = ba + 16*bi[i]; 1808e6580ceeSShri Abhyankar pj = bj + bi[i]; 1809e6580ceeSShri Abhyankar nz = bi[i+1] - bi[i]; 1810e6580ceeSShri Abhyankar 1811e6580ceeSShri Abhyankar /* Copy x block back into pv block */ 1812e6580ceeSShri Abhyankar for (j=0; j<nz; j++) { 1813e6580ceeSShri Abhyankar x = rtmp+16*((unsigned int)pj[j]); 1814e6580ceeSShri Abhyankar /* x = rtmp+4*pj[j]; */ 1815e6580ceeSShri Abhyankar 1816e6580ceeSShri Abhyankar SSE_INLINE_BEGIN_2(x,pv) 1817e6580ceeSShri Abhyankar /* Note: on future SSE architectures, STORE might be more efficient than STOREL/H */ 1818e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_0,XMM1) 1819e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_0,XMM1) 1820e6580ceeSShri Abhyankar 1821e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_2,XMM2) 1822e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_2,XMM2) 1823e6580ceeSShri Abhyankar 1824e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_4,XMM3) 1825e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_4,XMM3) 1826e6580ceeSShri Abhyankar 1827e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_6,XMM4) 1828e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_6,XMM4) 1829e6580ceeSShri Abhyankar 1830e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_8,XMM5) 1831e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_8,XMM5) 1832e6580ceeSShri Abhyankar 1833e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_10,XMM6) 1834e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_10,XMM6) 1835e6580ceeSShri Abhyankar 1836e6580ceeSShri Abhyankar SSE_LOADL_PS(SSE_ARG_1,FLOAT_12,XMM7) 1837e6580ceeSShri Abhyankar SSE_STOREL_PS(SSE_ARG_2,FLOAT_12,XMM7) 1838e6580ceeSShri Abhyankar 1839e6580ceeSShri Abhyankar SSE_LOADH_PS(SSE_ARG_1,FLOAT_14,XMM0) 1840e6580ceeSShri Abhyankar SSE_STOREH_PS(SSE_ARG_2,FLOAT_14,XMM0) 1841e6580ceeSShri Abhyankar SSE_INLINE_END_2; 1842e6580ceeSShri Abhyankar pv += 16; 1843e6580ceeSShri Abhyankar } 1844e6580ceeSShri Abhyankar /* invert diagonal block */ 1845e6580ceeSShri Abhyankar w = ba + 16*diag_offset[i]; 1846e6580ceeSShri Abhyankar if (pivotinblocks) { 1847*96b95a6bSBarry Smith ierr = PetscKernel_A_gets_inverse_A_4(w,shift);CHKERRQ(ierr); 1848e6580ceeSShri Abhyankar } else { 1849*96b95a6bSBarry Smith ierr = PetscKernel_A_gets_inverse_A_4_nopivot(w);CHKERRQ(ierr); 1850e6580ceeSShri Abhyankar } 1851*96b95a6bSBarry Smith /* ierr = PetscKernel_A_gets_inverse_A_4_SSE(w);CHKERRQ(ierr); */ 1852e6580ceeSShri Abhyankar /* Note: Using Kramer's rule, flop count below might be infairly high or low? */ 1853e6580ceeSShri Abhyankar } 1854e6580ceeSShri Abhyankar 1855e6580ceeSShri Abhyankar ierr = PetscFree(rtmp);CHKERRQ(ierr); 1856e6580ceeSShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering_SSE; 1857e6580ceeSShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_4_NaturalOrdering_SSE; 1858e6580ceeSShri Abhyankar C->assembled = PETSC_TRUE; 1859766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*bs*bs2*b->mbs);CHKERRQ(ierr); 1860e6580ceeSShri Abhyankar /* Flop Count from inverting diagonal blocks */ 1861e6580ceeSShri Abhyankar SSE_SCOPE_END; 1862e6580ceeSShri Abhyankar PetscFunctionReturn(0); 1863e6580ceeSShri Abhyankar } 1864e6580ceeSShri Abhyankar 1865e6580ceeSShri Abhyankar #endif 1866