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