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