xref: /petsc/src/mat/impls/baij/seq/dgefa4.c (revision e32f2f54e699d0aa6e733466c00da7e34666fe5e)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
34224c193SBarry Smith /*
48a36c062SBarry Smith        Inverts 4 by 4 matrix using partial pivoting.
571c5468dSBarry Smith 
671c5468dSBarry Smith        Used by the sparse factorization routines in
7dd882469SBarry Smith      src/mat/impls/baij/seq
871c5468dSBarry Smith 
971c5468dSBarry Smith        This is a combination of the Linpack routines
1071c5468dSBarry Smith     dgefa() and dgedi() specialized for a size of 4.
1171c5468dSBarry Smith 
124224c193SBarry Smith */
13d382aafbSBarry Smith #include "petscsys.h"
144224c193SBarry Smith 
154a2ae208SSatish Balay #undef __FUNCT__
164a2ae208SSatish Balay #define __FUNCT__ "Kernel_A_gets_inverse_A_4"
1762bba022SBarry Smith PetscErrorCode Kernel_A_gets_inverse_A_4(MatScalar *a,PetscReal shift)
184224c193SBarry Smith {
19690b6cddSBarry Smith     PetscInt   i__2,i__3,kp1,j,k,l,ll,i,ipvt[4],kb,k3;
20690b6cddSBarry Smith     PetscInt   k4,j3;
21b48ee343SBarry Smith     MatScalar  *aa,*ax,*ay,work[16],stmp;
22329f5518SBarry Smith     MatReal    tmp,max;
234224c193SBarry Smith 
244224c193SBarry Smith /*     gaussian elimination with partial pivoting */
254224c193SBarry Smith 
263a40ed3dSBarry Smith     PetscFunctionBegin;
271bfe5b3dSBarry Smith     shift = .25*shift*(PetscAbsScalar(a[0]) + PetscAbsScalar(a[5]) + PetscAbsScalar(a[10]) + PetscAbsScalar(a[15]));
284224c193SBarry Smith     /* Parameter adjustments */
298a36c062SBarry Smith     a       -= 5;
304224c193SBarry Smith 
318a36c062SBarry Smith     for (k = 1; k <= 3; ++k) {
324224c193SBarry Smith         kp1 = k + 1;
338a36c062SBarry Smith         k3  = 4*k;
344224c193SBarry Smith         k4  = k3 + k;
354224c193SBarry Smith /*        find l = pivot index */
364224c193SBarry Smith 
37ed33f8a5SSatish Balay         i__2 = 5 - k;
384224c193SBarry Smith         aa = &a[k4];
394224c193SBarry Smith         max = PetscAbsScalar(aa[0]);
404224c193SBarry Smith         l = 1;
414224c193SBarry Smith         for (ll=1; ll<i__2; ll++) {
424224c193SBarry Smith           tmp = PetscAbsScalar(aa[ll]);
434224c193SBarry Smith           if (tmp > max) { max = tmp; l = ll+1;}
444224c193SBarry Smith         }
454224c193SBarry Smith         l       += k - 1;
46da10e913SBarry Smith         ipvt[k-1] = l;
474224c193SBarry Smith 
485b8514ebSBarry Smith         if (a[l + k3] == 0.0) {
491bfe5b3dSBarry Smith           if (shift == 0.0) {
50*e32f2f54SBarry Smith 	    SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
511bfe5b3dSBarry Smith   	  } else {
521bfe5b3dSBarry Smith             /* SHIFT is applied to SINGLE diagonal entry; does this make any sense? */
531bfe5b3dSBarry Smith   	    a[l + k3] = shift;
541bfe5b3dSBarry Smith   	  }
554224c193SBarry Smith         }
564224c193SBarry Smith 
574224c193SBarry Smith /*           interchange if necessary */
584224c193SBarry Smith 
594224c193SBarry Smith         if (l != k) {
604224c193SBarry Smith           stmp      = a[l + k3];
614224c193SBarry Smith           a[l + k3] = a[k4];
624224c193SBarry Smith           a[k4]     = stmp;
634224c193SBarry Smith         }
644224c193SBarry Smith 
654224c193SBarry Smith /*           compute multipliers */
664224c193SBarry Smith 
674224c193SBarry Smith         stmp = -1. / a[k4];
688a36c062SBarry Smith         i__2 = 4 - k;
694224c193SBarry Smith         aa = &a[1 + k4];
704224c193SBarry Smith         for (ll=0; ll<i__2; ll++) {
714224c193SBarry Smith           aa[ll] *= stmp;
724224c193SBarry Smith         }
734224c193SBarry Smith 
744224c193SBarry Smith /*           row elimination with column indexing */
754224c193SBarry Smith 
764224c193SBarry Smith         ax = &a[k4+1];
778a36c062SBarry Smith         for (j = kp1; j <= 4; ++j) {
788a36c062SBarry Smith             j3   = 4*j;
794224c193SBarry Smith             stmp = a[l + j3];
804224c193SBarry Smith             if (l != k) {
814224c193SBarry Smith               a[l + j3] = a[k + j3];
824224c193SBarry Smith               a[k + j3] = stmp;
834224c193SBarry Smith             }
844224c193SBarry Smith 
858a36c062SBarry Smith             i__3 = 4 - k;
864224c193SBarry Smith             ay = &a[1+k+j3];
874224c193SBarry Smith             for (ll=0; ll<i__3; ll++) {
884224c193SBarry Smith               ay[ll] += stmp*ax[ll];
894224c193SBarry Smith             }
904224c193SBarry Smith         }
914224c193SBarry Smith     }
92da10e913SBarry Smith     ipvt[3] = 4;
935b8514ebSBarry Smith     if (a[20] == 0.0) {
94*e32f2f54SBarry Smith       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",3);
954224c193SBarry Smith     }
964224c193SBarry Smith 
974224c193SBarry Smith     /*
984224c193SBarry Smith          Now form the inverse
994224c193SBarry Smith     */
1004224c193SBarry Smith 
1014224c193SBarry Smith    /*     compute inverse(u) */
1024224c193SBarry Smith 
1038a36c062SBarry Smith     for (k = 1; k <= 4; ++k) {
1048a36c062SBarry Smith         k3    = 4*k;
1054224c193SBarry Smith         k4    = k3 + k;
1064224c193SBarry Smith         a[k4] = 1.0 / a[k4];
1074224c193SBarry Smith         stmp  = -a[k4];
1084224c193SBarry Smith         i__2  = k - 1;
1094224c193SBarry Smith         aa    = &a[k3 + 1];
1104224c193SBarry Smith         for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
1114224c193SBarry Smith         kp1 = k + 1;
1128a36c062SBarry Smith         if (4 < kp1) continue;
1134224c193SBarry Smith         ax = aa;
1148a36c062SBarry Smith         for (j = kp1; j <= 4; ++j) {
1158a36c062SBarry Smith             j3        = 4*j;
1164224c193SBarry Smith             stmp      = a[k + j3];
1174224c193SBarry Smith             a[k + j3] = 0.0;
1184224c193SBarry Smith             ay        = &a[j3 + 1];
1194224c193SBarry Smith             for (ll=0; ll<k; ll++) {
1204224c193SBarry Smith               ay[ll] += stmp*ax[ll];
1214224c193SBarry Smith             }
1224224c193SBarry Smith         }
1234224c193SBarry Smith     }
1244224c193SBarry Smith 
1254224c193SBarry Smith    /*    form inverse(u)*inverse(l) */
1264224c193SBarry Smith 
1278a36c062SBarry Smith     for (kb = 1; kb <= 3; ++kb) {
1288a36c062SBarry Smith         k   = 4 - kb;
1298a36c062SBarry Smith         k3  = 4*k;
1304224c193SBarry Smith         kp1 = k + 1;
1314224c193SBarry Smith         aa  = a + k3;
1328a36c062SBarry Smith         for (i = kp1; i <= 4; ++i) {
133b48ee343SBarry Smith             work[i-1] = aa[i];
1344224c193SBarry Smith             aa[i]   = 0.0;
1354224c193SBarry Smith         }
1368a36c062SBarry Smith         for (j = kp1; j <= 4; ++j) {
137b48ee343SBarry Smith             stmp  = work[j-1];
1388a36c062SBarry Smith             ax    = &a[4*j + 1];
1394224c193SBarry Smith             ay    = &a[k3 + 1];
1404224c193SBarry Smith             ay[0] += stmp*ax[0];
1414224c193SBarry Smith             ay[1] += stmp*ax[1];
1424224c193SBarry Smith             ay[2] += stmp*ax[2];
1438a36c062SBarry Smith             ay[3] += stmp*ax[3];
1444224c193SBarry Smith         }
145da10e913SBarry Smith         l = ipvt[k-1];
1464224c193SBarry Smith         if (l != k) {
1474224c193SBarry Smith             ax = &a[k3 + 1];
1488a36c062SBarry Smith             ay = &a[4*l + 1];
1494224c193SBarry Smith             stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
1504224c193SBarry Smith             stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
1514224c193SBarry Smith             stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
1528a36c062SBarry Smith             stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp;
1534224c193SBarry Smith         }
1544224c193SBarry Smith     }
1553a40ed3dSBarry Smith     PetscFunctionReturn(0);
1564224c193SBarry Smith }
1574224c193SBarry Smith 
1589b445006SKris Buschelman #if defined(PETSC_HAVE_SSE)
1599b445006SKris Buschelman #include PETSC_HAVE_SSE
16099148eceSKris Buschelman 
16199148eceSKris Buschelman #undef __FUNCT__
1627059b3dcSKris Buschelman #define __FUNCT__ "Kernel_A_gets_inverse_A_4_SSE"
163dfbe8321SBarry Smith PetscErrorCode Kernel_A_gets_inverse_A_4_SSE(float *a)
16499148eceSKris Buschelman {
16599148eceSKris Buschelman   /*
1669b445006SKris Buschelman      This routine is converted from Intel's Small Matrix Library.
16799148eceSKris Buschelman      See: Streaming SIMD Extensions -- Inverse of 4x4 Matrix
16899148eceSKris Buschelman      Order Number: 245043-001
16999148eceSKris Buschelman      March 1999
17099148eceSKris Buschelman      http://www.intel.com
17199148eceSKris Buschelman 
17299148eceSKris Buschelman      Inverse of a 4x4 matrix via Kramer's Rule:
17399148eceSKris Buschelman      bool Invert4x4(SMLXMatrix &);
17499148eceSKris Buschelman   */
17599148eceSKris Buschelman   PetscFunctionBegin;
1769b445006SKris Buschelman 
1779b445006SKris Buschelman   SSE_SCOPE_BEGIN;
1789b445006SKris Buschelman     SSE_INLINE_BEGIN_1(a)
1799b445006SKris Buschelman 
18099148eceSKris Buschelman /* ----------------------------------------------- */
1819b445006SKris Buschelman 
1829b445006SKris Buschelman       SSE_LOADL_PS(SSE_ARG_1,FLOAT_0,XMM0)
1839b445006SKris Buschelman       SSE_LOADH_PS(SSE_ARG_1,FLOAT_4,XMM0)
1849b445006SKris Buschelman 
1859b445006SKris Buschelman       SSE_LOADL_PS(SSE_ARG_1,FLOAT_8,XMM5)
1869b445006SKris Buschelman       SSE_LOADH_PS(SSE_ARG_1,FLOAT_12,XMM5)
1879b445006SKris Buschelman 
1889b445006SKris Buschelman       SSE_COPY_PS(XMM3,XMM0)
1899b445006SKris Buschelman       SSE_SHUFFLE(XMM3,XMM5,0x88)
1909b445006SKris Buschelman 
1919b445006SKris Buschelman       SSE_SHUFFLE(XMM5,XMM0,0xDD)
1929b445006SKris Buschelman 
1939b445006SKris Buschelman       SSE_LOADL_PS(SSE_ARG_1,FLOAT_2,XMM0)
1949b445006SKris Buschelman       SSE_LOADH_PS(SSE_ARG_1,FLOAT_6,XMM0)
1959b445006SKris Buschelman 
1969b445006SKris Buschelman       SSE_LOADL_PS(SSE_ARG_1,FLOAT_10,XMM6)
1979b445006SKris Buschelman       SSE_LOADH_PS(SSE_ARG_1,FLOAT_14,XMM6)
1989b445006SKris Buschelman 
1999b445006SKris Buschelman       SSE_COPY_PS(XMM4,XMM0)
2009b445006SKris Buschelman       SSE_SHUFFLE(XMM4,XMM6,0x88)
2019b445006SKris Buschelman 
2029b445006SKris Buschelman       SSE_SHUFFLE(XMM6,XMM0,0xDD)
2039b445006SKris Buschelman 
20499148eceSKris Buschelman /* ----------------------------------------------- */
2059b445006SKris Buschelman 
2069b445006SKris Buschelman       SSE_COPY_PS(XMM7,XMM4)
2079b445006SKris Buschelman       SSE_MULT_PS(XMM7,XMM6)
2089b445006SKris Buschelman 
2099b445006SKris Buschelman       SSE_SHUFFLE(XMM7,XMM7,0xB1)
2109b445006SKris Buschelman 
2119b445006SKris Buschelman       SSE_COPY_PS(XMM0,XMM5)
2129b445006SKris Buschelman       SSE_MULT_PS(XMM0,XMM7)
2139b445006SKris Buschelman 
2149b445006SKris Buschelman       SSE_COPY_PS(XMM2,XMM3)
2159b445006SKris Buschelman       SSE_MULT_PS(XMM2,XMM7)
2169b445006SKris Buschelman 
2179b445006SKris Buschelman       SSE_SHUFFLE(XMM7,XMM7,0x4E)
2189b445006SKris Buschelman 
2199b445006SKris Buschelman       SSE_COPY_PS(XMM1,XMM5)
2209b445006SKris Buschelman       SSE_MULT_PS(XMM1,XMM7)
2219b445006SKris Buschelman       SSE_SUB_PS(XMM1,XMM0)
2229b445006SKris Buschelman 
2239b445006SKris Buschelman       SSE_MULT_PS(XMM7,XMM3)
2249b445006SKris Buschelman       SSE_SUB_PS(XMM7,XMM2)
2259b445006SKris Buschelman 
2269b445006SKris Buschelman       SSE_SHUFFLE(XMM7,XMM7,0x4E)
2279b445006SKris Buschelman       SSE_STORE_PS(SSE_ARG_1,FLOAT_4,XMM7)
2289b445006SKris Buschelman 
22999148eceSKris Buschelman /* ----------------------------------------------- */
2309b445006SKris Buschelman 
2319b445006SKris Buschelman       SSE_COPY_PS(XMM0,XMM5)
2329b445006SKris Buschelman       SSE_MULT_PS(XMM0,XMM4)
2339b445006SKris Buschelman 
2349b445006SKris Buschelman       SSE_SHUFFLE(XMM0,XMM0,0xB1)
2359b445006SKris Buschelman 
2369b445006SKris Buschelman       SSE_COPY_PS(XMM2,XMM6)
2379b445006SKris Buschelman       SSE_MULT_PS(XMM2,XMM0)
2389b445006SKris Buschelman       SSE_ADD_PS(XMM2,XMM1)
2399b445006SKris Buschelman 
2409b445006SKris Buschelman       SSE_COPY_PS(XMM7,XMM3)
2419b445006SKris Buschelman       SSE_MULT_PS(XMM7,XMM0)
2429b445006SKris Buschelman 
2439b445006SKris Buschelman       SSE_SHUFFLE(XMM0,XMM0,0x4E)
2449b445006SKris Buschelman 
2459b445006SKris Buschelman       SSE_COPY_PS(XMM1,XMM6)
2469b445006SKris Buschelman       SSE_MULT_PS(XMM1,XMM0)
2479b445006SKris Buschelman       SSE_SUB_PS(XMM2,XMM1)
2489b445006SKris Buschelman 
2499b445006SKris Buschelman       SSE_MULT_PS(XMM0,XMM3)
2509b445006SKris Buschelman       SSE_SUB_PS(XMM0,XMM7)
2519b445006SKris Buschelman 
2529b445006SKris Buschelman       SSE_SHUFFLE(XMM0,XMM0,0x4E)
2539b445006SKris Buschelman       SSE_STORE_PS(SSE_ARG_1,FLOAT_12,XMM0)
2549b445006SKris Buschelman 
25599148eceSKris Buschelman       /* ----------------------------------------------- */
2569b445006SKris Buschelman 
2579b445006SKris Buschelman       SSE_COPY_PS(XMM7,XMM5)
2589b445006SKris Buschelman       SSE_SHUFFLE(XMM7,XMM5,0x4E)
2599b445006SKris Buschelman       SSE_MULT_PS(XMM7,XMM6)
2609b445006SKris Buschelman 
2619b445006SKris Buschelman       SSE_SHUFFLE(XMM7,XMM7,0xB1)
2629b445006SKris Buschelman 
2639b445006SKris Buschelman       SSE_SHUFFLE(XMM4,XMM4,0x4E)
2649b445006SKris Buschelman 
2659b445006SKris Buschelman       SSE_COPY_PS(XMM0,XMM4)
2669b445006SKris Buschelman       SSE_MULT_PS(XMM0,XMM7)
2679b445006SKris Buschelman       SSE_ADD_PS(XMM0,XMM2)
2689b445006SKris Buschelman 
2699b445006SKris Buschelman       SSE_COPY_PS(XMM2,XMM3)
2709b445006SKris Buschelman       SSE_MULT_PS(XMM2,XMM7)
2719b445006SKris Buschelman 
2729b445006SKris Buschelman       SSE_SHUFFLE(XMM7,XMM7,0x4E)
2739b445006SKris Buschelman 
2749b445006SKris Buschelman       SSE_COPY_PS(XMM1,XMM4)
2759b445006SKris Buschelman       SSE_MULT_PS(XMM1,XMM7)
2769b445006SKris Buschelman       SSE_SUB_PS(XMM0,XMM1)
2779b445006SKris Buschelman       SSE_STORE_PS(SSE_ARG_1,FLOAT_0,XMM0)
2789b445006SKris Buschelman 
2799b445006SKris Buschelman       SSE_MULT_PS(XMM7,XMM3)
2809b445006SKris Buschelman       SSE_SUB_PS(XMM7,XMM2)
2819b445006SKris Buschelman 
2829b445006SKris Buschelman       SSE_SHUFFLE(XMM7,XMM7,0x4E)
2839b445006SKris Buschelman 
28499148eceSKris Buschelman       /* ----------------------------------------------- */
2859b445006SKris Buschelman 
2869b445006SKris Buschelman       SSE_COPY_PS(XMM1,XMM3)
2879b445006SKris Buschelman       SSE_MULT_PS(XMM1,XMM5)
2889b445006SKris Buschelman 
2899b445006SKris Buschelman       SSE_SHUFFLE(XMM1,XMM1,0xB1)
2909b445006SKris Buschelman 
2919b445006SKris Buschelman       SSE_COPY_PS(XMM0,XMM6)
2929b445006SKris Buschelman       SSE_MULT_PS(XMM0,XMM1)
2939b445006SKris Buschelman       SSE_ADD_PS(XMM0,XMM7)
2949b445006SKris Buschelman 
2959b445006SKris Buschelman       SSE_COPY_PS(XMM2,XMM4)
2969b445006SKris Buschelman       SSE_MULT_PS(XMM2,XMM1)
2979b445006SKris Buschelman       SSE_SUB_PS_M(XMM2,SSE_ARG_1,FLOAT_12)
2989b445006SKris Buschelman 
2999b445006SKris Buschelman       SSE_SHUFFLE(XMM1,XMM1,0x4E)
3009b445006SKris Buschelman 
3019b445006SKris Buschelman       SSE_COPY_PS(XMM7,XMM6)
3029b445006SKris Buschelman       SSE_MULT_PS(XMM7,XMM1)
3039b445006SKris Buschelman       SSE_SUB_PS(XMM7,XMM0)
3049b445006SKris Buschelman 
3059b445006SKris Buschelman       SSE_MULT_PS(XMM1,XMM4)
3069b445006SKris Buschelman       SSE_SUB_PS(XMM2,XMM1)
3079b445006SKris Buschelman       SSE_STORE_PS(SSE_ARG_1,FLOAT_12,XMM2)
3089b445006SKris Buschelman 
30999148eceSKris Buschelman       /* ----------------------------------------------- */
3109b445006SKris Buschelman 
3119b445006SKris Buschelman       SSE_COPY_PS(XMM1,XMM3)
3129b445006SKris Buschelman       SSE_MULT_PS(XMM1,XMM6)
3139b445006SKris Buschelman 
3149b445006SKris Buschelman       SSE_SHUFFLE(XMM1,XMM1,0xB1)
3159b445006SKris Buschelman 
3169b445006SKris Buschelman       SSE_COPY_PS(XMM2,XMM4)
3179b445006SKris Buschelman       SSE_MULT_PS(XMM2,XMM1)
3189b445006SKris Buschelman       SSE_LOAD_PS(SSE_ARG_1,FLOAT_4,XMM0)
3199b445006SKris Buschelman       SSE_SUB_PS(XMM0,XMM2)
3209b445006SKris Buschelman 
3219b445006SKris Buschelman       SSE_COPY_PS(XMM2,XMM5)
3229b445006SKris Buschelman       SSE_MULT_PS(XMM2,XMM1)
3239b445006SKris Buschelman       SSE_ADD_PS(XMM2,XMM7)
3249b445006SKris Buschelman 
3259b445006SKris Buschelman       SSE_SHUFFLE(XMM1,XMM1,0x4E)
3269b445006SKris Buschelman 
3279b445006SKris Buschelman       SSE_COPY_PS(XMM7,XMM4)
3289b445006SKris Buschelman       SSE_MULT_PS(XMM7,XMM1)
3299b445006SKris Buschelman       SSE_ADD_PS(XMM7,XMM0)
3309b445006SKris Buschelman 
3319b445006SKris Buschelman       SSE_MULT_PS(XMM1,XMM5)
3329b445006SKris Buschelman       SSE_SUB_PS(XMM2,XMM1)
3339b445006SKris Buschelman 
33499148eceSKris Buschelman       /* ----------------------------------------------- */
3359b445006SKris Buschelman 
3369b445006SKris Buschelman       SSE_MULT_PS(XMM4,XMM3)
3379b445006SKris Buschelman 
3389b445006SKris Buschelman       SSE_SHUFFLE(XMM4,XMM4,0xB1)
3399b445006SKris Buschelman 
3409b445006SKris Buschelman       SSE_COPY_PS(XMM1,XMM6)
3419b445006SKris Buschelman       SSE_MULT_PS(XMM1,XMM4)
3429b445006SKris Buschelman       SSE_ADD_PS(XMM1,XMM7)
3439b445006SKris Buschelman 
3449b445006SKris Buschelman       SSE_COPY_PS(XMM0,XMM5)
3459b445006SKris Buschelman       SSE_MULT_PS(XMM0,XMM4)
3469b445006SKris Buschelman       SSE_LOAD_PS(SSE_ARG_1,FLOAT_12,XMM7)
3479b445006SKris Buschelman       SSE_SUB_PS(XMM7,XMM0)
3489b445006SKris Buschelman 
3499b445006SKris Buschelman       SSE_SHUFFLE(XMM4,XMM4,0x4E)
3509b445006SKris Buschelman 
3519b445006SKris Buschelman       SSE_MULT_PS(XMM6,XMM4)
3529b445006SKris Buschelman       SSE_SUB_PS(XMM1,XMM6)
3539b445006SKris Buschelman 
3549b445006SKris Buschelman       SSE_MULT_PS(XMM5,XMM4)
3559b445006SKris Buschelman       SSE_ADD_PS(XMM5,XMM7)
3569b445006SKris Buschelman 
3579b445006SKris Buschelman       /* ----------------------------------------------- */
3589b445006SKris Buschelman 
3599b445006SKris Buschelman       SSE_LOAD_PS(SSE_ARG_1,FLOAT_0,XMM0)
3609b445006SKris Buschelman       SSE_MULT_PS(XMM3,XMM0)
3619b445006SKris Buschelman 
3629b445006SKris Buschelman       SSE_COPY_PS(XMM4,XMM3)
3639b445006SKris Buschelman       SSE_SHUFFLE(XMM4,XMM3,0x4E)
3649b445006SKris Buschelman       SSE_ADD_PS(XMM4,XMM3)
3659b445006SKris Buschelman 
3669b445006SKris Buschelman       SSE_COPY_PS(XMM6,XMM4)
3679b445006SKris Buschelman       SSE_SHUFFLE(XMM6,XMM4,0xB1)
3689b445006SKris Buschelman       SSE_ADD_SS(XMM6,XMM4)
3699b445006SKris Buschelman 
3709b445006SKris Buschelman       SSE_COPY_PS(XMM3,XMM6)
3719b445006SKris Buschelman       SSE_RECIP_SS(XMM3,XMM6)
3729b445006SKris Buschelman       SSE_COPY_SS(XMM4,XMM3)
3739b445006SKris Buschelman       SSE_ADD_SS(XMM4,XMM3)
3749b445006SKris Buschelman       SSE_MULT_SS(XMM3,XMM3)
3759b445006SKris Buschelman       SSE_MULT_SS(XMM6,XMM3)
3769b445006SKris Buschelman       SSE_SUB_SS(XMM4,XMM6)
3779b445006SKris Buschelman 
3789b445006SKris Buschelman       SSE_SHUFFLE(XMM4,XMM4,0x00)
3799b445006SKris Buschelman 
3809b445006SKris Buschelman       SSE_MULT_PS(XMM0,XMM4)
3819b445006SKris Buschelman       SSE_STOREL_PS(SSE_ARG_1,FLOAT_0,XMM0)
3829b445006SKris Buschelman       SSE_STOREH_PS(SSE_ARG_1,FLOAT_2,XMM0)
3839b445006SKris Buschelman 
3849b445006SKris Buschelman       SSE_MULT_PS(XMM1,XMM4)
3859b445006SKris Buschelman       SSE_STOREL_PS(SSE_ARG_1,FLOAT_4,XMM1)
3869b445006SKris Buschelman       SSE_STOREH_PS(SSE_ARG_1,FLOAT_6,XMM1)
3879b445006SKris Buschelman 
3889b445006SKris Buschelman       SSE_MULT_PS(XMM2,XMM4)
3899b445006SKris Buschelman       SSE_STOREL_PS(SSE_ARG_1,FLOAT_8,XMM2)
3909b445006SKris Buschelman       SSE_STOREH_PS(SSE_ARG_1,FLOAT_10,XMM2)
3919b445006SKris Buschelman 
3929b445006SKris Buschelman       SSE_MULT_PS(XMM4,XMM5)
3939b445006SKris Buschelman       SSE_STOREL_PS(SSE_ARG_1,FLOAT_12,XMM4)
3949b445006SKris Buschelman       SSE_STOREH_PS(SSE_ARG_1,FLOAT_14,XMM4)
3959b445006SKris Buschelman 
3969b445006SKris Buschelman       /* ----------------------------------------------- */
3979b445006SKris Buschelman 
3989b445006SKris Buschelman       SSE_INLINE_END_1;
3999b445006SKris Buschelman   SSE_SCOPE_END;
4009b445006SKris Buschelman 
40199148eceSKris Buschelman   PetscFunctionReturn(0);
40299148eceSKris Buschelman }
40399148eceSKris Buschelman 
40499148eceSKris Buschelman #endif
40599148eceSKris Buschelman 
40699148eceSKris Buschelman 
407