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