xref: /petsc/src/mat/impls/baij/seq/dgefa2.c (revision 943c8ff524c3b3228b734140abd933ffca3f4d05)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
384643e36SBarry Smith /*
484643e36SBarry Smith      Inverts 2 by 2 matrix using partial pivoting.
584643e36SBarry Smith 
684643e36SBarry Smith        Used by the sparse factorization routines in
7dd882469SBarry Smith      src/mat/impls/baij/seq
884643e36SBarry Smith 
984643e36SBarry Smith 
1084643e36SBarry Smith        This is a combination of the Linpack routines
1184643e36SBarry Smith     dgefa() and dgedi() specialized for a size of 2.
1284643e36SBarry Smith 
1384643e36SBarry Smith */
14d382aafbSBarry Smith #include "petscsys.h"
1584643e36SBarry Smith 
164a2ae208SSatish Balay #undef __FUNCT__
174a2ae208SSatish Balay #define __FUNCT__ "Kernel_A_gets_inverse_A_2"
1862bba022SBarry Smith PetscErrorCode Kernel_A_gets_inverse_A_2(MatScalar *a,PetscReal shift)
1984643e36SBarry Smith {
20690b6cddSBarry Smith     PetscInt   i__2,i__3,kp1,j,k,l,ll,i,ipvt[2],k3;
21690b6cddSBarry Smith     PetscInt   k4,j3;
22b48ee343SBarry Smith     MatScalar  *aa,*ax,*ay,work[4],stmp;
23329f5518SBarry Smith     MatReal    tmp,max;
2484643e36SBarry Smith 
2584643e36SBarry Smith /*     gaussian elimination with partial pivoting */
2684643e36SBarry Smith 
2784643e36SBarry Smith     PetscFunctionBegin;
28*943c8ff5SShri Abhyankar     shift = .25*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[3]));
2984643e36SBarry Smith     /* Parameter adjustments */
3084643e36SBarry Smith     a       -= 3;
3184643e36SBarry Smith 
3284643e36SBarry Smith     /*for (k = 1; k <= 1; ++k) {*/
3384643e36SBarry Smith         k   = 1;
3484643e36SBarry Smith 	kp1 = k + 1;
3584643e36SBarry Smith         k3  = 2*k;
3684643e36SBarry Smith         k4  = k3 + k;
3784643e36SBarry Smith /*        find l = pivot index */
3884643e36SBarry Smith 
39ed33f8a5SSatish Balay 	i__2 = 3 - k;
4084643e36SBarry Smith         aa = &a[k4];
4184643e36SBarry Smith         max = PetscAbsScalar(aa[0]);
4284643e36SBarry Smith         l = 1;
4384643e36SBarry Smith         for (ll=1; ll<i__2; ll++) {
4484643e36SBarry Smith           tmp = PetscAbsScalar(aa[ll]);
4584643e36SBarry Smith           if (tmp > max) { max = tmp; l = ll+1;}
4684643e36SBarry Smith         }
4784643e36SBarry Smith         l       += k - 1;
48b48ee343SBarry Smith 	ipvt[k-1] = l;
4984643e36SBarry Smith 
50*943c8ff5SShri Abhyankar 	if (a[l + k3] == 0.0) {
51*943c8ff5SShri Abhyankar 	  if (shift == 0.0) {
52*943c8ff5SShri Abhyankar 	    SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
53*943c8ff5SShri Abhyankar 	  } else {
54*943c8ff5SShri Abhyankar 	    a[l + k3] = shift;
55*943c8ff5SShri Abhyankar 	  }
56*943c8ff5SShri Abhyankar 	}
5784643e36SBarry Smith 
5884643e36SBarry Smith /*           interchange if necessary */
5984643e36SBarry Smith 
6084643e36SBarry Smith 	if (l != k) {
6184643e36SBarry Smith 	  stmp      = a[l + k3];
6284643e36SBarry Smith 	  a[l + k3] = a[k4];
6384643e36SBarry Smith 	  a[k4]     = stmp;
6484643e36SBarry Smith         }
6584643e36SBarry Smith 
6684643e36SBarry Smith /*           compute multipliers */
6784643e36SBarry Smith 
6884643e36SBarry Smith 	stmp = -1. / a[k4];
6984643e36SBarry Smith 	i__2 = 2 - k;
7084643e36SBarry Smith         aa = &a[1 + k4];
7184643e36SBarry Smith         for (ll=0; ll<i__2; ll++) {
7284643e36SBarry Smith           aa[ll] *= stmp;
7384643e36SBarry Smith         }
7484643e36SBarry Smith 
7584643e36SBarry Smith /*           row elimination with column indexing */
7684643e36SBarry Smith 
7784643e36SBarry Smith 	ax = &a[k4+1];
7884643e36SBarry Smith         for (j = kp1; j <= 2; ++j) {
7984643e36SBarry Smith             j3   = 2*j;
8084643e36SBarry Smith 	    stmp = a[l + j3];
8184643e36SBarry Smith 	    if (l != k) {
8284643e36SBarry Smith 	      a[l + j3] = a[k + j3];
8384643e36SBarry Smith 	      a[k + j3] = stmp;
8484643e36SBarry Smith             }
8584643e36SBarry Smith 
8684643e36SBarry Smith 	    i__3 = 2 - k;
8784643e36SBarry Smith             ay = &a[1+k+j3];
8884643e36SBarry Smith             for (ll=0; ll<i__3; ll++) {
8984643e36SBarry Smith               ay[ll] += stmp*ax[ll];
9084643e36SBarry Smith             }
9184643e36SBarry Smith 	}
9284643e36SBarry Smith     /*}*/
93b48ee343SBarry Smith     ipvt[1] = 2;
9465e19b50SBarry Smith     if (a[6] == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",1);
9584643e36SBarry Smith 
9684643e36SBarry Smith     /*
9784643e36SBarry Smith          Now form the inverse
9884643e36SBarry Smith     */
9984643e36SBarry Smith 
10084643e36SBarry Smith    /*     compute inverse(u) */
10184643e36SBarry Smith 
10284643e36SBarry Smith     for (k = 1; k <= 2; ++k) {
10384643e36SBarry Smith         k3    = 2*k;
10484643e36SBarry Smith         k4    = k3 + k;
10584643e36SBarry Smith 	a[k4] = 1.0 / a[k4];
10684643e36SBarry Smith 	stmp  = -a[k4];
10784643e36SBarry Smith 	i__2  = k - 1;
10884643e36SBarry Smith         aa    = &a[k3 + 1];
10984643e36SBarry Smith         for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
11084643e36SBarry Smith 	kp1 = k + 1;
11184643e36SBarry Smith 	if (2 < kp1) continue;
11284643e36SBarry Smith         ax = aa;
11384643e36SBarry Smith         for (j = kp1; j <= 2; ++j) {
11484643e36SBarry Smith             j3        = 2*j;
11584643e36SBarry Smith 	    stmp      = a[k + j3];
11684643e36SBarry Smith 	    a[k + j3] = 0.0;
11784643e36SBarry Smith             ay        = &a[j3 + 1];
11884643e36SBarry Smith             for (ll=0; ll<k; ll++) {
11984643e36SBarry Smith               ay[ll] += stmp*ax[ll];
12084643e36SBarry Smith             }
12184643e36SBarry Smith 	}
12284643e36SBarry Smith     }
12384643e36SBarry Smith 
12484643e36SBarry Smith    /*    form inverse(u)*inverse(l) */
12584643e36SBarry Smith 
12684643e36SBarry Smith     /*for (kb = 1; kb <= 1; ++kb) {*/
12784643e36SBarry Smith 
12884643e36SBarry Smith 	k   = 1;
12984643e36SBarry Smith         k3  = 2*k;
13084643e36SBarry Smith 	kp1 = k + 1;
13184643e36SBarry Smith         aa  = a + k3;
13284643e36SBarry Smith 	for (i = kp1; i <= 2; ++i) {
133b48ee343SBarry Smith             work[i-1] = aa[i];
13484643e36SBarry Smith 	    aa[i]   = 0.0;
13584643e36SBarry Smith 	}
13684643e36SBarry Smith 	for (j = kp1; j <= 2; ++j) {
137b48ee343SBarry Smith 	    stmp  = work[j-1];
13884643e36SBarry Smith             ax    = &a[2*j + 1];
13984643e36SBarry Smith             ay    = &a[k3 + 1];
14084643e36SBarry Smith             ay[0] += stmp*ax[0];
14184643e36SBarry Smith             ay[1] += stmp*ax[1];
14284643e36SBarry Smith 	}
143b48ee343SBarry Smith 	l = ipvt[k-1];
14484643e36SBarry Smith 	if (l != k) {
14584643e36SBarry Smith             ax = &a[k3 + 1];
14684643e36SBarry Smith             ay = &a[2*l + 1];
14784643e36SBarry Smith             stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
14884643e36SBarry Smith             stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
14984643e36SBarry Smith 	}
15084643e36SBarry Smith 
15184643e36SBarry Smith     PetscFunctionReturn(0);
15284643e36SBarry Smith }
15384643e36SBarry Smith 
1543dfa136dSBarry Smith #undef __FUNCT__
1553dfa136dSBarry Smith #define __FUNCT__ "Kernel_A_gets_inverse_A_9"
15662bba022SBarry Smith PetscErrorCode Kernel_A_gets_inverse_A_9(MatScalar *a,PetscReal shift)
1573dfa136dSBarry Smith {
1583dfa136dSBarry Smith     PetscInt   i__2,i__3,kp1,j,k,l,ll,i,ipvt[9],kb,k3;
1593dfa136dSBarry Smith     PetscInt   k4,j3;
1603dfa136dSBarry Smith     MatScalar  *aa,*ax,*ay,work[81],stmp;
1613dfa136dSBarry Smith     MatReal    tmp,max;
1623dfa136dSBarry Smith 
1633dfa136dSBarry Smith /*     gaussian elimination with partial pivoting */
1643dfa136dSBarry Smith 
1653dfa136dSBarry Smith     PetscFunctionBegin;
1663dfa136dSBarry Smith     /* Parameter adjustments */
1673dfa136dSBarry Smith     a       -= 10;
1683dfa136dSBarry Smith 
1693dfa136dSBarry Smith     for (k = 1; k <= 8; ++k) {
1703dfa136dSBarry Smith 	kp1 = k + 1;
1713dfa136dSBarry Smith         k3  = 9*k;
1723dfa136dSBarry Smith         k4  = k3 + k;
1733dfa136dSBarry Smith /*        find l = pivot index */
1743dfa136dSBarry Smith 
175c38ccd74SBarry Smith 	i__2 = 10 - k;
1763dfa136dSBarry Smith         aa = &a[k4];
1773dfa136dSBarry Smith         max = PetscAbsScalar(aa[0]);
1783dfa136dSBarry Smith         l = 1;
1793dfa136dSBarry Smith         for (ll=1; ll<i__2; ll++) {
1803dfa136dSBarry Smith           tmp = PetscAbsScalar(aa[ll]);
1813dfa136dSBarry Smith           if (tmp > max) { max = tmp; l = ll+1;}
1823dfa136dSBarry Smith         }
1833dfa136dSBarry Smith         l       += k - 1;
1843dfa136dSBarry Smith 	ipvt[k-1] = l;
1853dfa136dSBarry Smith 
18665e19b50SBarry Smith 	if (a[l + k3] == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
1873dfa136dSBarry Smith 
1883dfa136dSBarry Smith /*           interchange if necessary */
1893dfa136dSBarry Smith 
1903dfa136dSBarry Smith 	if (l != k) {
1913dfa136dSBarry Smith 	  stmp      = a[l + k3];
1923dfa136dSBarry Smith 	  a[l + k3] = a[k4];
1933dfa136dSBarry Smith 	  a[k4]     = stmp;
1943dfa136dSBarry Smith         }
1953dfa136dSBarry Smith 
1963dfa136dSBarry Smith /*           compute multipliers */
1973dfa136dSBarry Smith 
1983dfa136dSBarry Smith 	stmp = -1. / a[k4];
1993dfa136dSBarry Smith 	i__2 = 9 - k;
2003dfa136dSBarry Smith         aa = &a[1 + k4];
2013dfa136dSBarry Smith         for (ll=0; ll<i__2; ll++) {
2023dfa136dSBarry Smith           aa[ll] *= stmp;
2033dfa136dSBarry Smith         }
2043dfa136dSBarry Smith 
2053dfa136dSBarry Smith /*           row elimination with column indexing */
2063dfa136dSBarry Smith 
2073dfa136dSBarry Smith 	ax = &a[k4+1];
2083dfa136dSBarry Smith         for (j = kp1; j <= 9; ++j) {
2093dfa136dSBarry Smith             j3   = 9*j;
2103dfa136dSBarry Smith 	    stmp = a[l + j3];
2113dfa136dSBarry Smith 	    if (l != k) {
2123dfa136dSBarry Smith 	      a[l + j3] = a[k + j3];
2133dfa136dSBarry Smith 	      a[k + j3] = stmp;
2143dfa136dSBarry Smith             }
2153dfa136dSBarry Smith 
2163dfa136dSBarry Smith 	    i__3 = 9 - k;
2173dfa136dSBarry Smith             ay = &a[1+k+j3];
2183dfa136dSBarry Smith             for (ll=0; ll<i__3; ll++) {
2193dfa136dSBarry Smith               ay[ll] += stmp*ax[ll];
2203dfa136dSBarry Smith             }
2213dfa136dSBarry Smith 	}
2223dfa136dSBarry Smith     }
2233dfa136dSBarry Smith     ipvt[8] = 9;
22465e19b50SBarry Smith     if (a[90] == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",6);
2253dfa136dSBarry Smith 
2263dfa136dSBarry Smith     /*
2273dfa136dSBarry Smith          Now form the inverse
2283dfa136dSBarry Smith     */
2293dfa136dSBarry Smith 
2303dfa136dSBarry Smith    /*     compute inverse(u) */
2313dfa136dSBarry Smith 
2323dfa136dSBarry Smith     for (k = 1; k <= 9; ++k) {
2333dfa136dSBarry Smith         k3    = 9*k;
2343dfa136dSBarry Smith         k4    = k3 + k;
2353dfa136dSBarry Smith 	a[k4] = 1.0 / a[k4];
2363dfa136dSBarry Smith 	stmp  = -a[k4];
2373dfa136dSBarry Smith 	i__2  = k - 1;
2383dfa136dSBarry Smith         aa    = &a[k3 + 1];
2393dfa136dSBarry Smith         for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
2403dfa136dSBarry Smith 	kp1 = k + 1;
2413dfa136dSBarry Smith 	if (9 < kp1) continue;
2423dfa136dSBarry Smith         ax = aa;
2433dfa136dSBarry Smith         for (j = kp1; j <= 9; ++j) {
2443dfa136dSBarry Smith             j3        = 9*j;
2453dfa136dSBarry Smith 	    stmp      = a[k + j3];
2463dfa136dSBarry Smith 	    a[k + j3] = 0.0;
2473dfa136dSBarry Smith             ay        = &a[j3 + 1];
2483dfa136dSBarry Smith             for (ll=0; ll<k; ll++) {
2493dfa136dSBarry Smith               ay[ll] += stmp*ax[ll];
2503dfa136dSBarry Smith             }
2513dfa136dSBarry Smith 	}
2523dfa136dSBarry Smith     }
2533dfa136dSBarry Smith 
2543dfa136dSBarry Smith    /*    form inverse(u)*inverse(l) */
2553dfa136dSBarry Smith 
2563dfa136dSBarry Smith     for (kb = 1; kb <= 8; ++kb) {
2573dfa136dSBarry Smith 	k   = 9 - kb;
2583dfa136dSBarry Smith         k3  = 9*k;
2593dfa136dSBarry Smith 	kp1 = k + 1;
2603dfa136dSBarry Smith         aa  = a + k3;
2613dfa136dSBarry Smith 	for (i = kp1; i <= 9; ++i) {
2623dfa136dSBarry Smith             work[i-1] = aa[i];
2633dfa136dSBarry Smith 	    aa[i]   = 0.0;
2643dfa136dSBarry Smith 	}
2653dfa136dSBarry Smith 	for (j = kp1; j <= 9; ++j) {
2663dfa136dSBarry Smith 	    stmp  = work[j-1];
2673dfa136dSBarry Smith             ax    = &a[9*j + 1];
2683dfa136dSBarry Smith             ay    = &a[k3 + 1];
2693dfa136dSBarry Smith             ay[0] += stmp*ax[0];
2703dfa136dSBarry Smith             ay[1] += stmp*ax[1];
2713dfa136dSBarry Smith             ay[2] += stmp*ax[2];
2723dfa136dSBarry Smith             ay[3] += stmp*ax[3];
2733dfa136dSBarry Smith             ay[4] += stmp*ax[4];
2743dfa136dSBarry Smith             ay[5] += stmp*ax[5];
2753dfa136dSBarry Smith             ay[6] += stmp*ax[6];
2763dfa136dSBarry Smith             ay[7] += stmp*ax[7];
2773dfa136dSBarry Smith             ay[8] += stmp*ax[8];
2783dfa136dSBarry Smith 	}
2793dfa136dSBarry Smith 	l = ipvt[k-1];
2803dfa136dSBarry Smith 	if (l != k) {
2813dfa136dSBarry Smith             ax = &a[k3 + 1];
2823dfa136dSBarry Smith             ay = &a[9*l + 1];
2833dfa136dSBarry Smith             stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
2843dfa136dSBarry Smith             stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
2853dfa136dSBarry Smith             stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
2863dfa136dSBarry Smith             stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp;
2873dfa136dSBarry Smith             stmp = ax[4]; ax[4] = ay[4]; ay[4] = stmp;
2883dfa136dSBarry Smith             stmp = ax[5]; ax[5] = ay[5]; ay[5] = stmp;
2893dfa136dSBarry Smith             stmp = ax[6]; ax[6] = ay[6]; ay[6] = stmp;
2903dfa136dSBarry Smith             stmp = ax[7]; ax[7] = ay[7]; ay[7] = stmp;
2913dfa136dSBarry Smith             stmp = ax[8]; ax[8] = ay[8]; ay[8] = stmp;
2923dfa136dSBarry Smith 	}
2933dfa136dSBarry Smith     }
2943dfa136dSBarry Smith     PetscFunctionReturn(0);
2953dfa136dSBarry Smith }
29684643e36SBarry Smith 
29729a97285SShri Abhyankar /*
29829a97285SShri Abhyankar       Inverts 15 by 15 matrix using partial pivoting.
29929a97285SShri Abhyankar 
30029a97285SShri Abhyankar        Used by the sparse factorization routines in
30129a97285SShri Abhyankar      src/mat/impls/baij/seq
30229a97285SShri Abhyankar 
30329a97285SShri Abhyankar        This is a combination of the Linpack routines
30429a97285SShri Abhyankar     dgefa() and dgedi() specialized for a size of 15.
30529a97285SShri Abhyankar 
30629a97285SShri Abhyankar */
30729a97285SShri Abhyankar 
30829a97285SShri Abhyankar #undef __FUNCT__
30929a97285SShri Abhyankar #define __FUNCT__ "Kernel_A_gets_inverse_A_15"
310766f9fbaSBarry Smith PetscErrorCode Kernel_A_gets_inverse_A_15(MatScalar *a,PetscInt *ipvt,MatScalar *work,PetscReal shift)
31129a97285SShri Abhyankar {
312766f9fbaSBarry Smith     PetscInt         i__2,i__3,kp1,j,k,l,ll,i,kb,k3;
31329a97285SShri Abhyankar     PetscInt         k4,j3;
314766f9fbaSBarry Smith     MatScalar        *aa,*ax,*ay,stmp;
31529a97285SShri Abhyankar     MatReal          tmp,max;
31629a97285SShri Abhyankar 
31729a97285SShri Abhyankar /*     gaussian elimination with partial pivoting */
31829a97285SShri Abhyankar 
31929a97285SShri Abhyankar     PetscFunctionBegin;
32029a97285SShri Abhyankar     /* Parameter adjustments */
32129a97285SShri Abhyankar     a       -= 16;
32229a97285SShri Abhyankar 
32329a97285SShri Abhyankar     for (k = 1; k <= 14; ++k) {
32429a97285SShri Abhyankar 	kp1 = k + 1;
32529a97285SShri Abhyankar         k3  = 15*k;
32629a97285SShri Abhyankar         k4  = k3 + k;
32729a97285SShri Abhyankar /*        find l = pivot index */
32829a97285SShri Abhyankar 
32929a97285SShri Abhyankar 	i__2 = 16 - k;
33029a97285SShri Abhyankar         aa = &a[k4];
33129a97285SShri Abhyankar         max = PetscAbsScalar(aa[0]);
33229a97285SShri Abhyankar         l = 1;
33329a97285SShri Abhyankar         for (ll=1; ll<i__2; ll++) {
33429a97285SShri Abhyankar           tmp = PetscAbsScalar(aa[ll]);
33529a97285SShri Abhyankar           if (tmp > max) { max = tmp; l = ll+1;}
33629a97285SShri Abhyankar         }
33729a97285SShri Abhyankar         l       += k - 1;
33829a97285SShri Abhyankar 	ipvt[k-1] = l;
33929a97285SShri Abhyankar 
34065e19b50SBarry Smith 	if (a[l + k3] == 0.0)  SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
34129a97285SShri Abhyankar 
34229a97285SShri Abhyankar /*           interchange if necessary */
34329a97285SShri Abhyankar 
34429a97285SShri Abhyankar 	if (l != k) {
34529a97285SShri Abhyankar 	  stmp      = a[l + k3];
34629a97285SShri Abhyankar 	  a[l + k3] = a[k4];
34729a97285SShri Abhyankar 	  a[k4]     = stmp;
34829a97285SShri Abhyankar         }
34929a97285SShri Abhyankar 
35029a97285SShri Abhyankar /*           compute multipliers */
35129a97285SShri Abhyankar 
35229a97285SShri Abhyankar 	stmp = -1. / a[k4];
35329a97285SShri Abhyankar 	i__2 = 15 - k;
35429a97285SShri Abhyankar         aa = &a[1 + k4];
35529a97285SShri Abhyankar         for (ll=0; ll<i__2; ll++) {
35629a97285SShri Abhyankar           aa[ll] *= stmp;
35729a97285SShri Abhyankar         }
35829a97285SShri Abhyankar 
35929a97285SShri Abhyankar /*           row elimination with column indexing */
36029a97285SShri Abhyankar 
36129a97285SShri Abhyankar 	ax = &a[k4+1];
36229a97285SShri Abhyankar         for (j = kp1; j <= 15; ++j) {
36329a97285SShri Abhyankar             j3   = 15*j;
36429a97285SShri Abhyankar 	    stmp = a[l + j3];
36529a97285SShri Abhyankar 	    if (l != k) {
36629a97285SShri Abhyankar 	      a[l + j3] = a[k + j3];
36729a97285SShri Abhyankar 	      a[k + j3] = stmp;
36829a97285SShri Abhyankar             }
36929a97285SShri Abhyankar 
37029a97285SShri Abhyankar 	    i__3 = 15 - k;
37129a97285SShri Abhyankar             ay = &a[1+k+j3];
37229a97285SShri Abhyankar             for (ll=0; ll<i__3; ll++) {
37329a97285SShri Abhyankar               ay[ll] += stmp*ax[ll];
37429a97285SShri Abhyankar             }
37529a97285SShri Abhyankar 	}
37629a97285SShri Abhyankar     }
37729a97285SShri Abhyankar     ipvt[14] = 15;
37865e19b50SBarry Smith     if (a[240] == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",6);
37929a97285SShri Abhyankar 
38029a97285SShri Abhyankar     /*
38129a97285SShri Abhyankar          Now form the inverse
38229a97285SShri Abhyankar     */
38329a97285SShri Abhyankar 
38429a97285SShri Abhyankar    /*     compute inverse(u) */
38529a97285SShri Abhyankar 
38629a97285SShri Abhyankar     for (k = 1; k <= 15; ++k) {
38729a97285SShri Abhyankar         k3    = 15*k;
38829a97285SShri Abhyankar         k4    = k3 + k;
38929a97285SShri Abhyankar 	a[k4] = 1.0 / a[k4];
39029a97285SShri Abhyankar 	stmp  = -a[k4];
39129a97285SShri Abhyankar 	i__2  = k - 1;
39229a97285SShri Abhyankar         aa    = &a[k3 + 1];
39329a97285SShri Abhyankar         for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
39429a97285SShri Abhyankar 	kp1 = k + 1;
39529a97285SShri Abhyankar 	if (15 < kp1) continue;
39629a97285SShri Abhyankar         ax = aa;
39729a97285SShri Abhyankar         for (j = kp1; j <= 15; ++j) {
39829a97285SShri Abhyankar             j3        = 15*j;
39929a97285SShri Abhyankar 	    stmp      = a[k + j3];
40029a97285SShri Abhyankar 	    a[k + j3] = 0.0;
40129a97285SShri Abhyankar             ay        = &a[j3 + 1];
40229a97285SShri Abhyankar             for (ll=0; ll<k; ll++) {
40329a97285SShri Abhyankar               ay[ll] += stmp*ax[ll];
40429a97285SShri Abhyankar             }
40529a97285SShri Abhyankar 	}
40629a97285SShri Abhyankar     }
40729a97285SShri Abhyankar 
40829a97285SShri Abhyankar    /*    form inverse(u)*inverse(l) */
40929a97285SShri Abhyankar 
41029a97285SShri Abhyankar     for (kb = 1; kb <= 14; ++kb) {
41129a97285SShri Abhyankar 	k   = 15 - kb;
41229a97285SShri Abhyankar         k3  = 15*k;
41329a97285SShri Abhyankar 	kp1 = k + 1;
41429a97285SShri Abhyankar         aa  = a + k3;
41529a97285SShri Abhyankar 	for (i = kp1; i <= 15; ++i) {
41629a97285SShri Abhyankar             work[i-1] = aa[i];
41729a97285SShri Abhyankar 	    aa[i]   = 0.0;
41829a97285SShri Abhyankar 	}
41929a97285SShri Abhyankar 	for (j = kp1; j <= 15; ++j) {
42029a97285SShri Abhyankar 	    stmp  = work[j-1];
42129a97285SShri Abhyankar             ax    = &a[15*j + 1];
42229a97285SShri Abhyankar             ay    = &a[k3 + 1];
42329a97285SShri Abhyankar             ay[0]  += stmp*ax[0];
42429a97285SShri Abhyankar             ay[1]  += stmp*ax[1];
42529a97285SShri Abhyankar             ay[2]  += stmp*ax[2];
42629a97285SShri Abhyankar             ay[3]  += stmp*ax[3];
42729a97285SShri Abhyankar             ay[4]  += stmp*ax[4];
42829a97285SShri Abhyankar             ay[5]  += stmp*ax[5];
42929a97285SShri Abhyankar             ay[6]  += stmp*ax[6];
43029a97285SShri Abhyankar 	    ay[7]  += stmp*ax[7];
43129a97285SShri Abhyankar             ay[8]  += stmp*ax[8];
43229a97285SShri Abhyankar             ay[9]  += stmp*ax[9];
43329a97285SShri Abhyankar             ay[10] += stmp*ax[10];
43429a97285SShri Abhyankar             ay[11] += stmp*ax[11];
43529a97285SShri Abhyankar             ay[12] += stmp*ax[12];
43629a97285SShri Abhyankar             ay[13] += stmp*ax[13];
43729a97285SShri Abhyankar 	    ay[14] += stmp*ax[14];
43829a97285SShri Abhyankar 	}
43929a97285SShri Abhyankar 	l = ipvt[k-1];
44029a97285SShri Abhyankar 	if (l != k) {
44129a97285SShri Abhyankar             ax = &a[k3 + 1];
44229a97285SShri Abhyankar             ay = &a[15*l + 1];
44329a97285SShri Abhyankar             stmp = ax[0];  ax[0]  = ay[0];  ay[0]  = stmp;
44429a97285SShri Abhyankar             stmp = ax[1];  ax[1]  = ay[1];  ay[1]  = stmp;
44529a97285SShri Abhyankar             stmp = ax[2];  ax[2]  = ay[2];  ay[2]  = stmp;
44629a97285SShri Abhyankar             stmp = ax[3];  ax[3]  = ay[3];  ay[3]  = stmp;
44729a97285SShri Abhyankar             stmp = ax[4];  ax[4]  = ay[4];  ay[4]  = stmp;
44829a97285SShri Abhyankar             stmp = ax[5];  ax[5]  = ay[5];  ay[5]  = stmp;
44929a97285SShri Abhyankar             stmp = ax[6];  ax[6]  = ay[6];  ay[6]  = stmp;
45029a97285SShri Abhyankar 	    stmp = ax[7];  ax[7]  = ay[7];  ay[7]  = stmp;
45129a97285SShri Abhyankar             stmp = ax[8];  ax[8]  = ay[8];  ay[8]  = stmp;
45229a97285SShri Abhyankar             stmp = ax[9];  ax[9]  = ay[9];  ay[9]  = stmp;
45329a97285SShri Abhyankar             stmp = ax[10]; ax[10] = ay[10]; ay[10] = stmp;
45429a97285SShri Abhyankar             stmp = ax[11]; ax[11] = ay[11]; ay[11] = stmp;
45529a97285SShri Abhyankar             stmp = ax[12]; ax[12] = ay[12]; ay[12] = stmp;
45629a97285SShri Abhyankar             stmp = ax[13]; ax[13] = ay[13]; ay[13] = stmp;
45729a97285SShri Abhyankar 	    stmp = ax[14]; ax[14] = ay[14]; ay[14] = stmp;
45829a97285SShri Abhyankar 	}
45929a97285SShri Abhyankar     }
46029a97285SShri Abhyankar     PetscFunctionReturn(0);
46129a97285SShri Abhyankar }
462