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