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