xref: /petsc/src/mat/impls/baij/seq/dgefa4.c (revision 26fbe8dc1a3c99fb8dddfa572c8c6b3b4ce3ca53)
1be1d678aSKris Buschelman 
24224c193SBarry Smith /*
38a36c062SBarry Smith        Inverts 4 by 4 matrix using partial pivoting.
471c5468dSBarry Smith 
571c5468dSBarry Smith        Used by the sparse factorization routines in
6dd882469SBarry Smith      src/mat/impls/baij/seq
771c5468dSBarry Smith 
871c5468dSBarry Smith        This is a combination of the Linpack routines
971c5468dSBarry Smith     dgefa() and dgedi() specialized for a size of 4.
1071c5468dSBarry Smith 
114224c193SBarry Smith */
12c6db04a5SJed Brown #include <petscsys.h>
134224c193SBarry Smith 
144a2ae208SSatish Balay #undef __FUNCT__
1596b95a6bSBarry Smith #define __FUNCT__ "PetscKernel_A_gets_inverse_A_4"
1696b95a6bSBarry Smith PetscErrorCode PetscKernel_A_gets_inverse_A_4(MatScalar *a,PetscReal shift)
174224c193SBarry Smith {
18690b6cddSBarry Smith   PetscInt  i__2,i__3,kp1,j,k,l,ll,i,ipvt[4],kb,k3;
19690b6cddSBarry Smith   PetscInt  k4,j3;
20b48ee343SBarry Smith   MatScalar *aa,*ax,*ay,work[16],stmp;
21329f5518SBarry Smith   MatReal   tmp,max;
224224c193SBarry Smith 
234224c193SBarry Smith /*     gaussian elimination with partial pivoting */
244224c193SBarry Smith 
253a40ed3dSBarry Smith   PetscFunctionBegin;
260daa89e9SHong Zhang   shift = .25*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[5]) + PetscAbsScalar(a[10]) + PetscAbsScalar(a[15]));
274224c193SBarry Smith   /* Parameter adjustments */
288a36c062SBarry Smith   a -= 5;
294224c193SBarry Smith 
308a36c062SBarry Smith   for (k = 1; k <= 3; ++k) {
314224c193SBarry Smith     kp1 = k + 1;
328a36c062SBarry Smith     k3  = 4*k;
334224c193SBarry Smith     k4  = k3 + k;
344224c193SBarry Smith /*        find l = pivot index */
354224c193SBarry Smith 
36ed33f8a5SSatish Balay     i__2 = 5 - k;
374224c193SBarry Smith     aa   = &a[k4];
384224c193SBarry Smith     max  = PetscAbsScalar(aa[0]);
394224c193SBarry Smith     l    = 1;
404224c193SBarry Smith     for (ll=1; ll<i__2; ll++) {
414224c193SBarry Smith       tmp = PetscAbsScalar(aa[ll]);
424224c193SBarry Smith       if (tmp > max) { max = tmp; l = ll+1;}
434224c193SBarry Smith     }
444224c193SBarry Smith     l        += k - 1;
45da10e913SBarry Smith     ipvt[k-1] = l;
464224c193SBarry Smith 
475b8514ebSBarry Smith     if (a[l + k3] == 0.0) {
48*26fbe8dcSKarl Rupp       if (shift == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
49*26fbe8dcSKarl Rupp       else {
501bfe5b3dSBarry Smith         /* SHIFT is applied to SINGLE diagonal entry; does this make any sense? */
511bfe5b3dSBarry Smith         a[l + k3] = shift;
521bfe5b3dSBarry Smith       }
534224c193SBarry Smith     }
544224c193SBarry Smith 
554224c193SBarry Smith /*           interchange if necessary */
564224c193SBarry Smith 
574224c193SBarry Smith     if (l != k) {
584224c193SBarry Smith       stmp      = a[l + k3];
594224c193SBarry Smith       a[l + k3] = a[k4];
604224c193SBarry Smith       a[k4]     = stmp;
614224c193SBarry Smith     }
624224c193SBarry Smith 
634224c193SBarry Smith /*           compute multipliers */
644224c193SBarry Smith 
654224c193SBarry Smith     stmp = -1. / a[k4];
668a36c062SBarry Smith     i__2 = 4 - k;
674224c193SBarry Smith     aa   = &a[1 + k4];
68*26fbe8dcSKarl Rupp     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
694224c193SBarry Smith 
704224c193SBarry Smith /*           row elimination with column indexing */
714224c193SBarry Smith 
724224c193SBarry Smith     ax = &a[k4+1];
738a36c062SBarry Smith     for (j = kp1; j <= 4; ++j) {
748a36c062SBarry Smith       j3   = 4*j;
754224c193SBarry Smith       stmp = a[l + j3];
764224c193SBarry Smith       if (l != k) {
774224c193SBarry Smith         a[l + j3] = a[k + j3];
784224c193SBarry Smith         a[k + j3] = stmp;
794224c193SBarry Smith       }
804224c193SBarry Smith 
818a36c062SBarry Smith       i__3 = 4 - k;
824224c193SBarry Smith       ay   = &a[1+k+j3];
83*26fbe8dcSKarl Rupp       for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll];
844224c193SBarry Smith     }
854224c193SBarry Smith   }
86da10e913SBarry Smith   ipvt[3] = 4;
8765e19b50SBarry Smith   if (a[20] == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",3);
884224c193SBarry Smith 
894224c193SBarry Smith   /*
904224c193SBarry Smith         Now form the inverse
914224c193SBarry Smith   */
924224c193SBarry Smith 
934224c193SBarry Smith   /*     compute inverse(u) */
944224c193SBarry Smith 
958a36c062SBarry Smith   for (k = 1; k <= 4; ++k) {
968a36c062SBarry Smith     k3    = 4*k;
974224c193SBarry Smith     k4    = k3 + k;
984224c193SBarry Smith     a[k4] = 1.0 / a[k4];
994224c193SBarry Smith     stmp  = -a[k4];
1004224c193SBarry Smith     i__2  = k - 1;
1014224c193SBarry Smith     aa    = &a[k3 + 1];
1024224c193SBarry Smith     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
1034224c193SBarry Smith     kp1 = k + 1;
1048a36c062SBarry Smith     if (4 < kp1) continue;
1054224c193SBarry Smith     ax = aa;
1068a36c062SBarry Smith     for (j = kp1; j <= 4; ++j) {
1078a36c062SBarry Smith       j3        = 4*j;
1084224c193SBarry Smith       stmp      = a[k + j3];
1094224c193SBarry Smith       a[k + j3] = 0.0;
1104224c193SBarry Smith       ay        = &a[j3 + 1];
111*26fbe8dcSKarl Rupp       for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll];
1124224c193SBarry Smith     }
1134224c193SBarry Smith   }
1144224c193SBarry Smith 
1154224c193SBarry Smith   /*    form inverse(u)*inverse(l) */
1164224c193SBarry Smith 
1178a36c062SBarry Smith   for (kb = 1; kb <= 3; ++kb) {
1188a36c062SBarry Smith     k   = 4 - kb;
1198a36c062SBarry Smith     k3  = 4*k;
1204224c193SBarry Smith     kp1 = k + 1;
1214224c193SBarry Smith     aa  = a + k3;
1228a36c062SBarry Smith     for (i = kp1; i <= 4; ++i) {
123b48ee343SBarry Smith       work[i-1] = aa[i];
1244224c193SBarry Smith       aa[i]     = 0.0;
1254224c193SBarry Smith     }
1268a36c062SBarry Smith     for (j = kp1; j <= 4; ++j) {
127b48ee343SBarry Smith       stmp   = work[j-1];
1288a36c062SBarry Smith       ax     = &a[4*j + 1];
1294224c193SBarry Smith       ay     = &a[k3 + 1];
1304224c193SBarry Smith       ay[0] += stmp*ax[0];
1314224c193SBarry Smith       ay[1] += stmp*ax[1];
1324224c193SBarry Smith       ay[2] += stmp*ax[2];
1338a36c062SBarry Smith       ay[3] += stmp*ax[3];
1344224c193SBarry Smith     }
135da10e913SBarry Smith     l = ipvt[k-1];
1364224c193SBarry Smith     if (l != k) {
1374224c193SBarry Smith       ax   = &a[k3 + 1];
1388a36c062SBarry Smith       ay   = &a[4*l + 1];
1394224c193SBarry Smith       stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
1404224c193SBarry Smith       stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
1414224c193SBarry Smith       stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
1428a36c062SBarry Smith       stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp;
1434224c193SBarry Smith     }
1444224c193SBarry Smith   }
1453a40ed3dSBarry Smith   PetscFunctionReturn(0);
1464224c193SBarry Smith }
1474224c193SBarry Smith 
1489b445006SKris Buschelman #if defined(PETSC_HAVE_SSE)
1499b445006SKris Buschelman #include PETSC_HAVE_SSE
15099148eceSKris Buschelman 
15199148eceSKris Buschelman #undef __FUNCT__
15296b95a6bSBarry Smith #define __FUNCT__ "PetscKernel_A_gets_inverse_A_4_SSE"
15396b95a6bSBarry Smith PetscErrorCode PetscKernel_A_gets_inverse_A_4_SSE(float *a)
15499148eceSKris Buschelman {
15599148eceSKris Buschelman   /*
1569b445006SKris Buschelman      This routine is converted from Intel's Small Matrix Library.
15799148eceSKris Buschelman      See: Streaming SIMD Extensions -- Inverse of 4x4 Matrix
15899148eceSKris Buschelman      Order Number: 245043-001
15999148eceSKris Buschelman      March 1999
16099148eceSKris Buschelman      http://www.intel.com
16199148eceSKris Buschelman 
16299148eceSKris Buschelman      Inverse of a 4x4 matrix via Kramer's Rule:
16399148eceSKris Buschelman      bool Invert4x4(SMLXMatrix &);
16499148eceSKris Buschelman   */
16599148eceSKris Buschelman   PetscFunctionBegin;
1669b445006SKris Buschelman   SSE_SCOPE_BEGIN;
1679b445006SKris Buschelman   SSE_INLINE_BEGIN_1(a)
1689b445006SKris Buschelman 
16999148eceSKris Buschelman /* ----------------------------------------------- */
1709b445006SKris Buschelman 
1719b445006SKris Buschelman   SSE_LOADL_PS(SSE_ARG_1,FLOAT_0,XMM0)
1729b445006SKris Buschelman   SSE_LOADH_PS(SSE_ARG_1,FLOAT_4,XMM0)
1739b445006SKris Buschelman 
1749b445006SKris Buschelman   SSE_LOADL_PS(SSE_ARG_1,FLOAT_8,XMM5)
1759b445006SKris Buschelman   SSE_LOADH_PS(SSE_ARG_1,FLOAT_12,XMM5)
1769b445006SKris Buschelman 
1779b445006SKris Buschelman   SSE_COPY_PS(XMM3,XMM0)
1789b445006SKris Buschelman   SSE_SHUFFLE(XMM3,XMM5,0x88)
1799b445006SKris Buschelman 
1809b445006SKris Buschelman   SSE_SHUFFLE(XMM5,XMM0,0xDD)
1819b445006SKris Buschelman 
1829b445006SKris Buschelman   SSE_LOADL_PS(SSE_ARG_1,FLOAT_2,XMM0)
1839b445006SKris Buschelman   SSE_LOADH_PS(SSE_ARG_1,FLOAT_6,XMM0)
1849b445006SKris Buschelman 
1859b445006SKris Buschelman   SSE_LOADL_PS(SSE_ARG_1,FLOAT_10,XMM6)
1869b445006SKris Buschelman   SSE_LOADH_PS(SSE_ARG_1,FLOAT_14,XMM6)
1879b445006SKris Buschelman 
1889b445006SKris Buschelman   SSE_COPY_PS(XMM4,XMM0)
1899b445006SKris Buschelman   SSE_SHUFFLE(XMM4,XMM6,0x88)
1909b445006SKris Buschelman 
1919b445006SKris Buschelman   SSE_SHUFFLE(XMM6,XMM0,0xDD)
1929b445006SKris Buschelman 
19399148eceSKris Buschelman /* ----------------------------------------------- */
1949b445006SKris Buschelman 
1959b445006SKris Buschelman   SSE_COPY_PS(XMM7,XMM4)
1969b445006SKris Buschelman   SSE_MULT_PS(XMM7,XMM6)
1979b445006SKris Buschelman 
1989b445006SKris Buschelman   SSE_SHUFFLE(XMM7,XMM7,0xB1)
1999b445006SKris Buschelman 
2009b445006SKris Buschelman   SSE_COPY_PS(XMM0,XMM5)
2019b445006SKris Buschelman   SSE_MULT_PS(XMM0,XMM7)
2029b445006SKris Buschelman 
2039b445006SKris Buschelman   SSE_COPY_PS(XMM2,XMM3)
2049b445006SKris Buschelman   SSE_MULT_PS(XMM2,XMM7)
2059b445006SKris Buschelman 
2069b445006SKris Buschelman   SSE_SHUFFLE(XMM7,XMM7,0x4E)
2079b445006SKris Buschelman 
2089b445006SKris Buschelman   SSE_COPY_PS(XMM1,XMM5)
2099b445006SKris Buschelman   SSE_MULT_PS(XMM1,XMM7)
2109b445006SKris Buschelman   SSE_SUB_PS(XMM1,XMM0)
2119b445006SKris Buschelman 
2129b445006SKris Buschelman   SSE_MULT_PS(XMM7,XMM3)
2139b445006SKris Buschelman   SSE_SUB_PS(XMM7,XMM2)
2149b445006SKris Buschelman 
2159b445006SKris Buschelman   SSE_SHUFFLE(XMM7,XMM7,0x4E)
2169b445006SKris Buschelman   SSE_STORE_PS(SSE_ARG_1,FLOAT_4,XMM7)
2179b445006SKris Buschelman 
21899148eceSKris Buschelman /* ----------------------------------------------- */
2199b445006SKris Buschelman 
2209b445006SKris Buschelman   SSE_COPY_PS(XMM0,XMM5)
2219b445006SKris Buschelman   SSE_MULT_PS(XMM0,XMM4)
2229b445006SKris Buschelman 
2239b445006SKris Buschelman   SSE_SHUFFLE(XMM0,XMM0,0xB1)
2249b445006SKris Buschelman 
2259b445006SKris Buschelman   SSE_COPY_PS(XMM2,XMM6)
2269b445006SKris Buschelman   SSE_MULT_PS(XMM2,XMM0)
2279b445006SKris Buschelman   SSE_ADD_PS(XMM2,XMM1)
2289b445006SKris Buschelman 
2299b445006SKris Buschelman   SSE_COPY_PS(XMM7,XMM3)
2309b445006SKris Buschelman   SSE_MULT_PS(XMM7,XMM0)
2319b445006SKris Buschelman 
2329b445006SKris Buschelman   SSE_SHUFFLE(XMM0,XMM0,0x4E)
2339b445006SKris Buschelman 
2349b445006SKris Buschelman   SSE_COPY_PS(XMM1,XMM6)
2359b445006SKris Buschelman   SSE_MULT_PS(XMM1,XMM0)
2369b445006SKris Buschelman   SSE_SUB_PS(XMM2,XMM1)
2379b445006SKris Buschelman 
2389b445006SKris Buschelman   SSE_MULT_PS(XMM0,XMM3)
2399b445006SKris Buschelman   SSE_SUB_PS(XMM0,XMM7)
2409b445006SKris Buschelman 
2419b445006SKris Buschelman   SSE_SHUFFLE(XMM0,XMM0,0x4E)
2429b445006SKris Buschelman   SSE_STORE_PS(SSE_ARG_1,FLOAT_12,XMM0)
2439b445006SKris Buschelman 
24499148eceSKris Buschelman   /* ----------------------------------------------- */
2459b445006SKris Buschelman 
2469b445006SKris Buschelman   SSE_COPY_PS(XMM7,XMM5)
2479b445006SKris Buschelman   SSE_SHUFFLE(XMM7,XMM5,0x4E)
2489b445006SKris Buschelman   SSE_MULT_PS(XMM7,XMM6)
2499b445006SKris Buschelman 
2509b445006SKris Buschelman   SSE_SHUFFLE(XMM7,XMM7,0xB1)
2519b445006SKris Buschelman 
2529b445006SKris Buschelman   SSE_SHUFFLE(XMM4,XMM4,0x4E)
2539b445006SKris Buschelman 
2549b445006SKris Buschelman   SSE_COPY_PS(XMM0,XMM4)
2559b445006SKris Buschelman   SSE_MULT_PS(XMM0,XMM7)
2569b445006SKris Buschelman   SSE_ADD_PS(XMM0,XMM2)
2579b445006SKris Buschelman 
2589b445006SKris Buschelman   SSE_COPY_PS(XMM2,XMM3)
2599b445006SKris Buschelman   SSE_MULT_PS(XMM2,XMM7)
2609b445006SKris Buschelman 
2619b445006SKris Buschelman   SSE_SHUFFLE(XMM7,XMM7,0x4E)
2629b445006SKris Buschelman 
2639b445006SKris Buschelman   SSE_COPY_PS(XMM1,XMM4)
2649b445006SKris Buschelman   SSE_MULT_PS(XMM1,XMM7)
2659b445006SKris Buschelman   SSE_SUB_PS(XMM0,XMM1)
2669b445006SKris Buschelman   SSE_STORE_PS(SSE_ARG_1,FLOAT_0,XMM0)
2679b445006SKris Buschelman 
2689b445006SKris Buschelman   SSE_MULT_PS(XMM7,XMM3)
2699b445006SKris Buschelman   SSE_SUB_PS(XMM7,XMM2)
2709b445006SKris Buschelman 
2719b445006SKris Buschelman   SSE_SHUFFLE(XMM7,XMM7,0x4E)
2729b445006SKris Buschelman 
27399148eceSKris Buschelman   /* ----------------------------------------------- */
2749b445006SKris Buschelman 
2759b445006SKris Buschelman   SSE_COPY_PS(XMM1,XMM3)
2769b445006SKris Buschelman   SSE_MULT_PS(XMM1,XMM5)
2779b445006SKris Buschelman 
2789b445006SKris Buschelman   SSE_SHUFFLE(XMM1,XMM1,0xB1)
2799b445006SKris Buschelman 
2809b445006SKris Buschelman   SSE_COPY_PS(XMM0,XMM6)
2819b445006SKris Buschelman   SSE_MULT_PS(XMM0,XMM1)
2829b445006SKris Buschelman   SSE_ADD_PS(XMM0,XMM7)
2839b445006SKris Buschelman 
2849b445006SKris Buschelman   SSE_COPY_PS(XMM2,XMM4)
2859b445006SKris Buschelman   SSE_MULT_PS(XMM2,XMM1)
2869b445006SKris Buschelman   SSE_SUB_PS_M(XMM2,SSE_ARG_1,FLOAT_12)
2879b445006SKris Buschelman 
2889b445006SKris Buschelman   SSE_SHUFFLE(XMM1,XMM1,0x4E)
2899b445006SKris Buschelman 
2909b445006SKris Buschelman   SSE_COPY_PS(XMM7,XMM6)
2919b445006SKris Buschelman   SSE_MULT_PS(XMM7,XMM1)
2929b445006SKris Buschelman   SSE_SUB_PS(XMM7,XMM0)
2939b445006SKris Buschelman 
2949b445006SKris Buschelman   SSE_MULT_PS(XMM1,XMM4)
2959b445006SKris Buschelman   SSE_SUB_PS(XMM2,XMM1)
2969b445006SKris Buschelman   SSE_STORE_PS(SSE_ARG_1,FLOAT_12,XMM2)
2979b445006SKris Buschelman 
29899148eceSKris Buschelman   /* ----------------------------------------------- */
2999b445006SKris Buschelman 
3009b445006SKris Buschelman   SSE_COPY_PS(XMM1,XMM3)
3019b445006SKris Buschelman   SSE_MULT_PS(XMM1,XMM6)
3029b445006SKris Buschelman 
3039b445006SKris Buschelman   SSE_SHUFFLE(XMM1,XMM1,0xB1)
3049b445006SKris Buschelman 
3059b445006SKris Buschelman   SSE_COPY_PS(XMM2,XMM4)
3069b445006SKris Buschelman   SSE_MULT_PS(XMM2,XMM1)
3079b445006SKris Buschelman   SSE_LOAD_PS(SSE_ARG_1,FLOAT_4,XMM0)
3089b445006SKris Buschelman   SSE_SUB_PS(XMM0,XMM2)
3099b445006SKris Buschelman 
3109b445006SKris Buschelman   SSE_COPY_PS(XMM2,XMM5)
3119b445006SKris Buschelman   SSE_MULT_PS(XMM2,XMM1)
3129b445006SKris Buschelman   SSE_ADD_PS(XMM2,XMM7)
3139b445006SKris Buschelman 
3149b445006SKris Buschelman   SSE_SHUFFLE(XMM1,XMM1,0x4E)
3159b445006SKris Buschelman 
3169b445006SKris Buschelman   SSE_COPY_PS(XMM7,XMM4)
3179b445006SKris Buschelman   SSE_MULT_PS(XMM7,XMM1)
3189b445006SKris Buschelman   SSE_ADD_PS(XMM7,XMM0)
3199b445006SKris Buschelman 
3209b445006SKris Buschelman   SSE_MULT_PS(XMM1,XMM5)
3219b445006SKris Buschelman   SSE_SUB_PS(XMM2,XMM1)
3229b445006SKris Buschelman 
32399148eceSKris Buschelman   /* ----------------------------------------------- */
3249b445006SKris Buschelman 
3259b445006SKris Buschelman   SSE_MULT_PS(XMM4,XMM3)
3269b445006SKris Buschelman 
3279b445006SKris Buschelman   SSE_SHUFFLE(XMM4,XMM4,0xB1)
3289b445006SKris Buschelman 
3299b445006SKris Buschelman   SSE_COPY_PS(XMM1,XMM6)
3309b445006SKris Buschelman   SSE_MULT_PS(XMM1,XMM4)
3319b445006SKris Buschelman   SSE_ADD_PS(XMM1,XMM7)
3329b445006SKris Buschelman 
3339b445006SKris Buschelman   SSE_COPY_PS(XMM0,XMM5)
3349b445006SKris Buschelman   SSE_MULT_PS(XMM0,XMM4)
3359b445006SKris Buschelman   SSE_LOAD_PS(SSE_ARG_1,FLOAT_12,XMM7)
3369b445006SKris Buschelman   SSE_SUB_PS(XMM7,XMM0)
3379b445006SKris Buschelman 
3389b445006SKris Buschelman   SSE_SHUFFLE(XMM4,XMM4,0x4E)
3399b445006SKris Buschelman 
3409b445006SKris Buschelman   SSE_MULT_PS(XMM6,XMM4)
3419b445006SKris Buschelman   SSE_SUB_PS(XMM1,XMM6)
3429b445006SKris Buschelman 
3439b445006SKris Buschelman   SSE_MULT_PS(XMM5,XMM4)
3449b445006SKris Buschelman   SSE_ADD_PS(XMM5,XMM7)
3459b445006SKris Buschelman 
3469b445006SKris Buschelman   /* ----------------------------------------------- */
3479b445006SKris Buschelman 
3489b445006SKris Buschelman   SSE_LOAD_PS(SSE_ARG_1,FLOAT_0,XMM0)
3499b445006SKris Buschelman   SSE_MULT_PS(XMM3,XMM0)
3509b445006SKris Buschelman 
3519b445006SKris Buschelman   SSE_COPY_PS(XMM4,XMM3)
3529b445006SKris Buschelman   SSE_SHUFFLE(XMM4,XMM3,0x4E)
3539b445006SKris Buschelman   SSE_ADD_PS(XMM4,XMM3)
3549b445006SKris Buschelman 
3559b445006SKris Buschelman   SSE_COPY_PS(XMM6,XMM4)
3569b445006SKris Buschelman   SSE_SHUFFLE(XMM6,XMM4,0xB1)
3579b445006SKris Buschelman   SSE_ADD_SS(XMM6,XMM4)
3589b445006SKris Buschelman 
3599b445006SKris Buschelman   SSE_COPY_PS(XMM3,XMM6)
3609b445006SKris Buschelman   SSE_RECIP_SS(XMM3,XMM6)
3619b445006SKris Buschelman   SSE_COPY_SS(XMM4,XMM3)
3629b445006SKris Buschelman   SSE_ADD_SS(XMM4,XMM3)
3639b445006SKris Buschelman   SSE_MULT_SS(XMM3,XMM3)
3649b445006SKris Buschelman   SSE_MULT_SS(XMM6,XMM3)
3659b445006SKris Buschelman   SSE_SUB_SS(XMM4,XMM6)
3669b445006SKris Buschelman 
3679b445006SKris Buschelman   SSE_SHUFFLE(XMM4,XMM4,0x00)
3689b445006SKris Buschelman 
3699b445006SKris Buschelman   SSE_MULT_PS(XMM0,XMM4)
3709b445006SKris Buschelman   SSE_STOREL_PS(SSE_ARG_1,FLOAT_0,XMM0)
3719b445006SKris Buschelman   SSE_STOREH_PS(SSE_ARG_1,FLOAT_2,XMM0)
3729b445006SKris Buschelman 
3739b445006SKris Buschelman   SSE_MULT_PS(XMM1,XMM4)
3749b445006SKris Buschelman   SSE_STOREL_PS(SSE_ARG_1,FLOAT_4,XMM1)
3759b445006SKris Buschelman   SSE_STOREH_PS(SSE_ARG_1,FLOAT_6,XMM1)
3769b445006SKris Buschelman 
3779b445006SKris Buschelman   SSE_MULT_PS(XMM2,XMM4)
3789b445006SKris Buschelman   SSE_STOREL_PS(SSE_ARG_1,FLOAT_8,XMM2)
3799b445006SKris Buschelman   SSE_STOREH_PS(SSE_ARG_1,FLOAT_10,XMM2)
3809b445006SKris Buschelman 
3819b445006SKris Buschelman   SSE_MULT_PS(XMM4,XMM5)
3829b445006SKris Buschelman   SSE_STOREL_PS(SSE_ARG_1,FLOAT_12,XMM4)
3839b445006SKris Buschelman   SSE_STOREH_PS(SSE_ARG_1,FLOAT_14,XMM4)
3849b445006SKris Buschelman 
3859b445006SKris Buschelman   /* ----------------------------------------------- */
3869b445006SKris Buschelman 
3879b445006SKris Buschelman   SSE_INLINE_END_1;
3889b445006SKris Buschelman   SSE_SCOPE_END;
38999148eceSKris Buschelman   PetscFunctionReturn(0);
39099148eceSKris Buschelman }
39199148eceSKris Buschelman 
39299148eceSKris Buschelman #endif
39399148eceSKris Buschelman 
39499148eceSKris Buschelman 
395