xref: /petsc/src/mat/impls/baij/seq/dgefa2.c (revision 70d8d27c0c8a665b79be2909c201b24b4c33566d)
1be1d678aSKris Buschelman 
284643e36SBarry Smith /*
3c9b7c560SHong Zhang      Inverts 2 by 2 matrix using gaussian elimination with partial pivoting.
484643e36SBarry Smith 
584643e36SBarry Smith        Used by the sparse factorization routines in
6dd882469SBarry Smith      src/mat/impls/baij/seq
784643e36SBarry Smith 
884643e36SBarry Smith 
984643e36SBarry Smith        This is a combination of the Linpack routines
1084643e36SBarry Smith     dgefa() and dgedi() specialized for a size of 2.
1184643e36SBarry Smith 
1284643e36SBarry Smith */
13c6db04a5SJed Brown #include <petscsys.h>
1484643e36SBarry Smith 
154a2ae208SSatish Balay #undef __FUNCT__
1696b95a6bSBarry Smith #define __FUNCT__ "PetscKernel_A_gets_inverse_A_2"
172e92ee13SHong Zhang PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_2(MatScalar *a,PetscReal shift,PetscBool allowzeropivot,PetscBool *zeropivotdetected)
1884643e36SBarry Smith {
19690b6cddSBarry Smith   PetscInt  i__2,i__3,kp1,j,k,l,ll,i,ipvt[2],k3;
20690b6cddSBarry Smith   PetscInt  k4,j3;
21b48ee343SBarry Smith   MatScalar *aa,*ax,*ay,work[4],stmp;
22329f5518SBarry Smith   MatReal   tmp,max;
2384643e36SBarry Smith 
2484643e36SBarry Smith   PetscFunctionBegin;
25c80103daSHong Zhang   if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
26943c8ff5SShri Abhyankar   shift = .25*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[3]));
27c9b7c560SHong Zhang 
2884643e36SBarry Smith   /* Parameter adjustments */
2984643e36SBarry Smith   a -= 3;
3084643e36SBarry Smith 
3184643e36SBarry Smith   k   = 1;
3284643e36SBarry Smith   kp1 = k + 1;
3384643e36SBarry Smith   k3  = 2*k;
3484643e36SBarry Smith   k4  = k3 + k;
3584643e36SBarry Smith 
36c9b7c560SHong Zhang   /* find l = pivot index */
37ed33f8a5SSatish Balay   i__2 = 3 - k;
3884643e36SBarry Smith   aa   = &a[k4];
3984643e36SBarry Smith   max  = PetscAbsScalar(aa[0]);
4084643e36SBarry Smith   l    = 1;
4184643e36SBarry Smith   for (ll=1; ll<i__2; ll++) {
4284643e36SBarry Smith     tmp = PetscAbsScalar(aa[ll]);
4384643e36SBarry Smith     if (tmp > max) { max = tmp; l = ll+1;}
4484643e36SBarry Smith   }
4584643e36SBarry Smith   l        += k - 1;
46b48ee343SBarry Smith   ipvt[k-1] = l;
4784643e36SBarry Smith 
48943c8ff5SShri Abhyankar   if (a[l + k3] == 0.0) {
49ec1892c8SHong Zhang     if (shift == 0.0) {
50ec1892c8SHong Zhang       if (allowzeropivot) {
51ec1892c8SHong Zhang         PetscErrorCode ierr;
52ec1892c8SHong Zhang         ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",k-1);CHKERRQ(ierr);
53*70d8d27cSBarry Smith         if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
54ec1892c8SHong Zhang       } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
55ec1892c8SHong Zhang     } else {
56943c8ff5SShri Abhyankar       a[l + k3] = shift;
57943c8ff5SShri Abhyankar     }
58943c8ff5SShri Abhyankar   }
5984643e36SBarry Smith 
6084643e36SBarry Smith   /* interchange if necessary */
6184643e36SBarry Smith   if (l != k) {
6284643e36SBarry Smith     stmp      = a[l + k3];
6384643e36SBarry Smith     a[l + k3] = a[k4];
6484643e36SBarry Smith     a[k4]     = stmp;
6584643e36SBarry Smith   }
6684643e36SBarry Smith 
6784643e36SBarry Smith   /* compute multipliers */
6884643e36SBarry Smith   stmp = -1. / a[k4];
6984643e36SBarry Smith   i__2 = 2 - k;
7084643e36SBarry Smith   aa = &a[1 + k4];
7126fbe8dcSKarl Rupp   for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
7284643e36SBarry Smith 
7384643e36SBarry Smith   /* row elimination with column indexing */
7484643e36SBarry Smith   ax = &a[k4+1];
7584643e36SBarry Smith   for (j = kp1; j <= 2; ++j) {
7684643e36SBarry Smith     j3   = 2*j;
7784643e36SBarry Smith     stmp = a[l + j3];
7884643e36SBarry Smith     if (l != k) {
7984643e36SBarry Smith       a[l + j3] = a[k + j3];
8084643e36SBarry Smith       a[k + j3] = stmp;
8184643e36SBarry Smith     }
8284643e36SBarry Smith 
8384643e36SBarry Smith     i__3 = 2 - k;
8484643e36SBarry Smith     ay   = &a[1+k+j3];
8526fbe8dcSKarl Rupp     for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll];
8684643e36SBarry Smith   }
8726fbe8dcSKarl Rupp 
88b48ee343SBarry Smith   ipvt[1] = 2;
892e92ee13SHong Zhang   if (a[6] == 0.0) {
902e92ee13SHong Zhang     if (allowzeropivot) {
91ec1892c8SHong Zhang       PetscErrorCode ierr;
922e92ee13SHong Zhang       ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",1);CHKERRQ(ierr);
93*70d8d27cSBarry Smith       if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
942e92ee13SHong Zhang     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",1);
952e92ee13SHong Zhang   }
9684643e36SBarry Smith 
97c9b7c560SHong Zhang   /* Now form the inverse */
9884643e36SBarry Smith   /* compute inverse(u) */
9984643e36SBarry Smith   for (k = 1; k <= 2; ++k) {
10084643e36SBarry Smith     k3    = 2*k;
10184643e36SBarry Smith     k4    = k3 + k;
10284643e36SBarry Smith     a[k4] = 1.0 / a[k4];
10384643e36SBarry Smith     stmp  = -a[k4];
10484643e36SBarry Smith     i__2  = k - 1;
10584643e36SBarry Smith     aa    = &a[k3 + 1];
10684643e36SBarry Smith     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
10784643e36SBarry Smith     kp1 = k + 1;
10884643e36SBarry Smith     if (2 < kp1) continue;
10984643e36SBarry Smith     ax = aa;
11084643e36SBarry Smith     for (j = kp1; j <= 2; ++j) {
11184643e36SBarry Smith       j3        = 2*j;
11284643e36SBarry Smith       stmp      = a[k + j3];
11384643e36SBarry Smith       a[k + j3] = 0.0;
11484643e36SBarry Smith       ay        = &a[j3 + 1];
11526fbe8dcSKarl Rupp       for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll];
11684643e36SBarry Smith     }
11784643e36SBarry Smith   }
11884643e36SBarry Smith 
11984643e36SBarry Smith   /* form inverse(u)*inverse(l) */
12084643e36SBarry Smith   k   = 1;
12184643e36SBarry Smith   k3  = 2*k;
12284643e36SBarry Smith   kp1 = k + 1;
12384643e36SBarry Smith   aa  = a + k3;
12484643e36SBarry Smith   for (i = kp1; i <= 2; ++i) {
125b48ee343SBarry Smith     work[i-1] = aa[i];
12684643e36SBarry Smith     aa[i]     = 0.0;
12784643e36SBarry Smith   }
12884643e36SBarry Smith   for (j = kp1; j <= 2; ++j) {
129b48ee343SBarry Smith     stmp   = work[j-1];
13084643e36SBarry Smith     ax     = &a[2*j + 1];
13184643e36SBarry Smith     ay     = &a[k3 + 1];
13284643e36SBarry Smith     ay[0] += stmp*ax[0];
13384643e36SBarry Smith     ay[1] += stmp*ax[1];
13484643e36SBarry Smith   }
135b48ee343SBarry Smith   l = ipvt[k-1];
13684643e36SBarry Smith   if (l != k) {
13784643e36SBarry Smith     ax   = &a[k3 + 1];
13884643e36SBarry Smith     ay   = &a[2*l + 1];
13984643e36SBarry Smith     stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
14084643e36SBarry Smith     stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
14184643e36SBarry Smith   }
14284643e36SBarry Smith   PetscFunctionReturn(0);
14384643e36SBarry Smith }
14484643e36SBarry Smith 
145ec1892c8SHong Zhang /* gaussian elimination with partial pivoting */
1463dfa136dSBarry Smith #undef __FUNCT__
14796b95a6bSBarry Smith #define __FUNCT__ "PetscKernel_A_gets_inverse_A_9"
1482e92ee13SHong Zhang PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_9(MatScalar *a,PetscReal shift,PetscBool allowzeropivot,PetscBool *zeropivotdetected)
1493dfa136dSBarry Smith {
1503dfa136dSBarry Smith   PetscInt  i__2,i__3,kp1,j,k,l,ll,i,ipvt[9],kb,k3;
1513dfa136dSBarry Smith   PetscInt  k4,j3;
1523dfa136dSBarry Smith   MatScalar *aa,*ax,*ay,work[81],stmp;
1533dfa136dSBarry Smith   MatReal   tmp,max;
1543dfa136dSBarry Smith 
1553dfa136dSBarry Smith   PetscFunctionBegin;
156c80103daSHong Zhang   if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
157c80103daSHong Zhang 
1583dfa136dSBarry Smith   /* Parameter adjustments */
1593dfa136dSBarry Smith   a -= 10;
1603dfa136dSBarry Smith 
1613dfa136dSBarry Smith   for (k = 1; k <= 8; ++k) {
1623dfa136dSBarry Smith     kp1 = k + 1;
1633dfa136dSBarry Smith     k3  = 9*k;
1643dfa136dSBarry Smith     k4  = k3 + k;
1653dfa136dSBarry Smith 
166ec1892c8SHong Zhang     /* find l = pivot index */
167c38ccd74SBarry Smith     i__2 = 10 - k;
1683dfa136dSBarry Smith     aa   = &a[k4];
1693dfa136dSBarry Smith     max  = PetscAbsScalar(aa[0]);
1703dfa136dSBarry Smith     l    = 1;
1713dfa136dSBarry Smith     for (ll=1; ll<i__2; ll++) {
1723dfa136dSBarry Smith       tmp = PetscAbsScalar(aa[ll]);
1733dfa136dSBarry Smith       if (tmp > max) { max = tmp; l = ll+1;}
1743dfa136dSBarry Smith     }
1753dfa136dSBarry Smith     l        += k - 1;
1763dfa136dSBarry Smith     ipvt[k-1] = l;
1773dfa136dSBarry Smith 
178ec1892c8SHong Zhang     if (a[l + k3] == 0.0) {
179ec1892c8SHong Zhang       if (shift == 0.0) {
180ec1892c8SHong Zhang         if (allowzeropivot) {
181ec1892c8SHong Zhang           PetscErrorCode ierr;
182ec1892c8SHong Zhang           ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",k-1);CHKERRQ(ierr);
183*70d8d27cSBarry Smith           if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
184ec1892c8SHong Zhang         } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
185ec1892c8SHong Zhang       } else {
186ec1892c8SHong Zhang         a[l + k3] = shift;
187ec1892c8SHong Zhang       }
188ec1892c8SHong Zhang     }
1893dfa136dSBarry Smith 
1903dfa136dSBarry Smith     /* interchange if necessary */
1913dfa136dSBarry Smith     if (l != k) {
1923dfa136dSBarry Smith       stmp      = a[l + k3];
1933dfa136dSBarry Smith       a[l + k3] = a[k4];
1943dfa136dSBarry Smith       a[k4]     = stmp;
1953dfa136dSBarry Smith     }
1963dfa136dSBarry Smith 
1973dfa136dSBarry Smith     /* compute multipliers */
1983dfa136dSBarry Smith     stmp = -1. / a[k4];
1993dfa136dSBarry Smith     i__2 = 9 - k;
2003dfa136dSBarry Smith     aa = &a[1 + k4];
20126fbe8dcSKarl Rupp     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
2023dfa136dSBarry Smith 
2033dfa136dSBarry Smith     /* row elimination with column indexing */
2043dfa136dSBarry Smith     ax = &a[k4+1];
2053dfa136dSBarry Smith     for (j = kp1; j <= 9; ++j) {
2063dfa136dSBarry Smith       j3   = 9*j;
2073dfa136dSBarry Smith       stmp = a[l + j3];
2083dfa136dSBarry Smith       if (l != k) {
2093dfa136dSBarry Smith         a[l + j3] = a[k + j3];
2103dfa136dSBarry Smith         a[k + j3] = stmp;
2113dfa136dSBarry Smith       }
2123dfa136dSBarry Smith 
2133dfa136dSBarry Smith       i__3 = 9 - k;
2143dfa136dSBarry Smith       ay = &a[1+k+j3];
21526fbe8dcSKarl Rupp       for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll];
2163dfa136dSBarry Smith     }
2173dfa136dSBarry Smith   }
2183dfa136dSBarry Smith   ipvt[8] = 9;
2192e92ee13SHong Zhang   if (a[90] == 0.0) {
2202e92ee13SHong Zhang     if (allowzeropivot) {
221ec1892c8SHong Zhang       PetscErrorCode ierr;
2222e92ee13SHong Zhang       ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",6);CHKERRQ(ierr);
223*70d8d27cSBarry Smith       if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
2242e92ee13SHong Zhang     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",6);
2252e92ee13SHong Zhang   }
2263dfa136dSBarry Smith 
227ec1892c8SHong Zhang   /* Now form the inverse */
2283dfa136dSBarry Smith   /* compute inverse(u) */
2293dfa136dSBarry Smith   for (k = 1; k <= 9; ++k) {
2303dfa136dSBarry Smith     k3    = 9*k;
2313dfa136dSBarry Smith     k4    = k3 + k;
2323dfa136dSBarry Smith     a[k4] = 1.0 / a[k4];
2333dfa136dSBarry Smith     stmp  = -a[k4];
2343dfa136dSBarry Smith     i__2  = k - 1;
2353dfa136dSBarry Smith     aa    = &a[k3 + 1];
2363dfa136dSBarry Smith     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
2373dfa136dSBarry Smith     kp1 = k + 1;
2383dfa136dSBarry Smith     if (9 < kp1) continue;
2393dfa136dSBarry Smith     ax = aa;
2403dfa136dSBarry Smith     for (j = kp1; j <= 9; ++j) {
2413dfa136dSBarry Smith       j3        = 9*j;
2423dfa136dSBarry Smith       stmp      = a[k + j3];
2433dfa136dSBarry Smith       a[k + j3] = 0.0;
2443dfa136dSBarry Smith       ay        = &a[j3 + 1];
24526fbe8dcSKarl Rupp       for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll];
2463dfa136dSBarry Smith     }
2473dfa136dSBarry Smith   }
2483dfa136dSBarry Smith 
2493dfa136dSBarry Smith   /* form inverse(u)*inverse(l) */
2503dfa136dSBarry Smith   for (kb = 1; kb <= 8; ++kb) {
2513dfa136dSBarry Smith     k   = 9 - kb;
2523dfa136dSBarry Smith     k3  = 9*k;
2533dfa136dSBarry Smith     kp1 = k + 1;
2543dfa136dSBarry Smith     aa  = a + k3;
2553dfa136dSBarry Smith     for (i = kp1; i <= 9; ++i) {
2563dfa136dSBarry Smith       work[i-1] = aa[i];
2573dfa136dSBarry Smith       aa[i]     = 0.0;
2583dfa136dSBarry Smith     }
2593dfa136dSBarry Smith     for (j = kp1; j <= 9; ++j) {
2603dfa136dSBarry Smith       stmp   = work[j-1];
2613dfa136dSBarry Smith       ax     = &a[9*j + 1];
2623dfa136dSBarry Smith       ay     = &a[k3 + 1];
2633dfa136dSBarry Smith       ay[0] += stmp*ax[0];
2643dfa136dSBarry Smith       ay[1] += stmp*ax[1];
2653dfa136dSBarry Smith       ay[2] += stmp*ax[2];
2663dfa136dSBarry Smith       ay[3] += stmp*ax[3];
2673dfa136dSBarry Smith       ay[4] += stmp*ax[4];
2683dfa136dSBarry Smith       ay[5] += stmp*ax[5];
2693dfa136dSBarry Smith       ay[6] += stmp*ax[6];
2703dfa136dSBarry Smith       ay[7] += stmp*ax[7];
2713dfa136dSBarry Smith       ay[8] += stmp*ax[8];
2723dfa136dSBarry Smith     }
2733dfa136dSBarry Smith     l = ipvt[k-1];
2743dfa136dSBarry Smith     if (l != k) {
2753dfa136dSBarry Smith       ax   = &a[k3 + 1];
2763dfa136dSBarry Smith       ay   = &a[9*l + 1];
2773dfa136dSBarry Smith       stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
2783dfa136dSBarry Smith       stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
2793dfa136dSBarry Smith       stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
2803dfa136dSBarry Smith       stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp;
2813dfa136dSBarry Smith       stmp = ax[4]; ax[4] = ay[4]; ay[4] = stmp;
2823dfa136dSBarry Smith       stmp = ax[5]; ax[5] = ay[5]; ay[5] = stmp;
2833dfa136dSBarry Smith       stmp = ax[6]; ax[6] = ay[6]; ay[6] = stmp;
2843dfa136dSBarry Smith       stmp = ax[7]; ax[7] = ay[7]; ay[7] = stmp;
2853dfa136dSBarry Smith       stmp = ax[8]; ax[8] = ay[8]; ay[8] = stmp;
2863dfa136dSBarry Smith     }
2873dfa136dSBarry Smith   }
2883dfa136dSBarry Smith   PetscFunctionReturn(0);
2893dfa136dSBarry Smith }
29084643e36SBarry Smith 
29129a97285SShri Abhyankar /*
292ec1892c8SHong Zhang       Inverts 15 by 15 matrix using gaussian elimination with partial pivoting.
29329a97285SShri Abhyankar 
29429a97285SShri Abhyankar        Used by the sparse factorization routines in
29529a97285SShri Abhyankar      src/mat/impls/baij/seq
29629a97285SShri Abhyankar 
29729a97285SShri Abhyankar        This is a combination of the Linpack routines
29829a97285SShri Abhyankar     dgefa() and dgedi() specialized for a size of 15.
29929a97285SShri Abhyankar 
30029a97285SShri Abhyankar */
30129a97285SShri Abhyankar 
30229a97285SShri Abhyankar #undef __FUNCT__
30396b95a6bSBarry Smith #define __FUNCT__ "PetscKernel_A_gets_inverse_A_15"
304922032d7SHong Zhang PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_15(MatScalar *a,PetscInt *ipvt,MatScalar *work,PetscReal shift,PetscBool allowzeropivot,PetscBool *zeropivotdetected)
30529a97285SShri Abhyankar {
306766f9fbaSBarry Smith   PetscInt  i__2,i__3,kp1,j,k,l,ll,i,kb,k3;
30729a97285SShri Abhyankar   PetscInt  k4,j3;
308766f9fbaSBarry Smith   MatScalar *aa,*ax,*ay,stmp;
30929a97285SShri Abhyankar   MatReal   tmp,max;
31029a97285SShri Abhyankar 
31129a97285SShri Abhyankar   PetscFunctionBegin;
312c80103daSHong Zhang   if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
313c80103daSHong Zhang 
31429a97285SShri Abhyankar   /* Parameter adjustments */
31529a97285SShri Abhyankar   a -= 16;
31629a97285SShri Abhyankar 
31729a97285SShri Abhyankar   for (k = 1; k <= 14; ++k) {
31829a97285SShri Abhyankar     kp1 = k + 1;
31929a97285SShri Abhyankar     k3  = 15*k;
32029a97285SShri Abhyankar     k4  = k3 + k;
32129a97285SShri Abhyankar 
322ec1892c8SHong Zhang     /* find l = pivot index */
32329a97285SShri Abhyankar     i__2 = 16 - k;
32429a97285SShri Abhyankar     aa   = &a[k4];
32529a97285SShri Abhyankar     max  = PetscAbsScalar(aa[0]);
32629a97285SShri Abhyankar     l    = 1;
32729a97285SShri Abhyankar     for (ll=1; ll<i__2; ll++) {
32829a97285SShri Abhyankar       tmp = PetscAbsScalar(aa[ll]);
32929a97285SShri Abhyankar       if (tmp > max) { max = tmp; l = ll+1;}
33029a97285SShri Abhyankar     }
33129a97285SShri Abhyankar     l        += k - 1;
33229a97285SShri Abhyankar     ipvt[k-1] = l;
33329a97285SShri Abhyankar 
334ec1892c8SHong Zhang     if (a[l + k3] == 0.0) {
335ec1892c8SHong Zhang       if (shift == 0.0) {
336ec1892c8SHong Zhang         if (allowzeropivot) {
337ec1892c8SHong Zhang           PetscErrorCode ierr;
338ec1892c8SHong Zhang           ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",k-1);CHKERRQ(ierr);
339*70d8d27cSBarry Smith           if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
340ec1892c8SHong Zhang         } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
341ec1892c8SHong Zhang       } else {
342ec1892c8SHong Zhang         a[l + k3] = shift;
343ec1892c8SHong Zhang       }
344ec1892c8SHong Zhang     }
34529a97285SShri Abhyankar 
34629a97285SShri Abhyankar     /* interchange if necessary */
34729a97285SShri Abhyankar     if (l != k) {
34829a97285SShri Abhyankar       stmp      = a[l + k3];
34929a97285SShri Abhyankar       a[l + k3] = a[k4];
35029a97285SShri Abhyankar       a[k4]     = stmp;
35129a97285SShri Abhyankar     }
35229a97285SShri Abhyankar 
35329a97285SShri Abhyankar     /* compute multipliers */
35429a97285SShri Abhyankar     stmp = -1. / a[k4];
35529a97285SShri Abhyankar     i__2 = 15 - k;
35629a97285SShri Abhyankar     aa = &a[1 + k4];
35726fbe8dcSKarl Rupp     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
35829a97285SShri Abhyankar 
35929a97285SShri Abhyankar     /* row elimination with column indexing */
36029a97285SShri Abhyankar     ax = &a[k4+1];
36129a97285SShri Abhyankar     for (j = kp1; j <= 15; ++j) {
36229a97285SShri Abhyankar       j3   = 15*j;
36329a97285SShri Abhyankar       stmp = a[l + j3];
36429a97285SShri Abhyankar       if (l != k) {
36529a97285SShri Abhyankar         a[l + j3] = a[k + j3];
36629a97285SShri Abhyankar         a[k + j3] = stmp;
36729a97285SShri Abhyankar       }
36829a97285SShri Abhyankar 
36929a97285SShri Abhyankar       i__3 = 15 - k;
37029a97285SShri Abhyankar       ay = &a[1+k+j3];
37126fbe8dcSKarl Rupp       for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll];
37229a97285SShri Abhyankar     }
37329a97285SShri Abhyankar   }
37429a97285SShri Abhyankar   ipvt[14] = 15;
375922032d7SHong Zhang   if (a[240] == 0.0) {
376922032d7SHong Zhang     if (allowzeropivot) {
377ec1892c8SHong Zhang       PetscErrorCode ierr;
378922032d7SHong Zhang       ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",6);CHKERRQ(ierr);
379*70d8d27cSBarry Smith       if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
380922032d7SHong Zhang     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",6);
381922032d7SHong Zhang   }
38229a97285SShri Abhyankar 
383ec1892c8SHong Zhang   /* Now form the inverse */
38429a97285SShri Abhyankar   /* compute inverse(u) */
38529a97285SShri Abhyankar   for (k = 1; k <= 15; ++k) {
38629a97285SShri Abhyankar     k3    = 15*k;
38729a97285SShri Abhyankar     k4    = k3 + k;
38829a97285SShri Abhyankar     a[k4] = 1.0 / a[k4];
38929a97285SShri Abhyankar     stmp  = -a[k4];
39029a97285SShri Abhyankar     i__2  = k - 1;
39129a97285SShri Abhyankar     aa    = &a[k3 + 1];
39229a97285SShri Abhyankar     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
39329a97285SShri Abhyankar     kp1 = k + 1;
39429a97285SShri Abhyankar     if (15 < kp1) continue;
39529a97285SShri Abhyankar     ax = aa;
39629a97285SShri Abhyankar     for (j = kp1; j <= 15; ++j) {
39729a97285SShri Abhyankar       j3        = 15*j;
39829a97285SShri Abhyankar       stmp      = a[k + j3];
39929a97285SShri Abhyankar       a[k + j3] = 0.0;
40029a97285SShri Abhyankar       ay        = &a[j3 + 1];
40126fbe8dcSKarl Rupp       for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll];
40229a97285SShri Abhyankar     }
40329a97285SShri Abhyankar   }
40429a97285SShri Abhyankar 
40529a97285SShri Abhyankar   /* form inverse(u)*inverse(l) */
40629a97285SShri Abhyankar   for (kb = 1; kb <= 14; ++kb) {
40729a97285SShri Abhyankar     k   = 15 - kb;
40829a97285SShri Abhyankar     k3  = 15*k;
40929a97285SShri Abhyankar     kp1 = k + 1;
41029a97285SShri Abhyankar     aa  = a + k3;
41129a97285SShri Abhyankar     for (i = kp1; i <= 15; ++i) {
41229a97285SShri Abhyankar       work[i-1] = aa[i];
41329a97285SShri Abhyankar       aa[i]     = 0.0;
41429a97285SShri Abhyankar     }
41529a97285SShri Abhyankar     for (j = kp1; j <= 15; ++j) {
41629a97285SShri Abhyankar       stmp    = work[j-1];
41729a97285SShri Abhyankar       ax      = &a[15*j + 1];
41829a97285SShri Abhyankar       ay      = &a[k3 + 1];
41929a97285SShri Abhyankar       ay[0]  += stmp*ax[0];
42029a97285SShri Abhyankar       ay[1]  += stmp*ax[1];
42129a97285SShri Abhyankar       ay[2]  += stmp*ax[2];
42229a97285SShri Abhyankar       ay[3]  += stmp*ax[3];
42329a97285SShri Abhyankar       ay[4]  += stmp*ax[4];
42429a97285SShri Abhyankar       ay[5]  += stmp*ax[5];
42529a97285SShri Abhyankar       ay[6]  += stmp*ax[6];
42629a97285SShri Abhyankar       ay[7]  += stmp*ax[7];
42729a97285SShri Abhyankar       ay[8]  += stmp*ax[8];
42829a97285SShri Abhyankar       ay[9]  += stmp*ax[9];
42929a97285SShri Abhyankar       ay[10] += stmp*ax[10];
43029a97285SShri Abhyankar       ay[11] += stmp*ax[11];
43129a97285SShri Abhyankar       ay[12] += stmp*ax[12];
43229a97285SShri Abhyankar       ay[13] += stmp*ax[13];
43329a97285SShri Abhyankar       ay[14] += stmp*ax[14];
43429a97285SShri Abhyankar     }
43529a97285SShri Abhyankar     l = ipvt[k-1];
43629a97285SShri Abhyankar     if (l != k) {
43729a97285SShri Abhyankar       ax   = &a[k3 + 1];
43829a97285SShri Abhyankar       ay   = &a[15*l + 1];
43929a97285SShri Abhyankar       stmp = ax[0];  ax[0]  = ay[0];  ay[0]  = stmp;
44029a97285SShri Abhyankar       stmp = ax[1];  ax[1]  = ay[1];  ay[1]  = stmp;
44129a97285SShri Abhyankar       stmp = ax[2];  ax[2]  = ay[2];  ay[2]  = stmp;
44229a97285SShri Abhyankar       stmp = ax[3];  ax[3]  = ay[3];  ay[3]  = stmp;
44329a97285SShri Abhyankar       stmp = ax[4];  ax[4]  = ay[4];  ay[4]  = stmp;
44429a97285SShri Abhyankar       stmp = ax[5];  ax[5]  = ay[5];  ay[5]  = stmp;
44529a97285SShri Abhyankar       stmp = ax[6];  ax[6]  = ay[6];  ay[6]  = stmp;
44629a97285SShri Abhyankar       stmp = ax[7];  ax[7]  = ay[7];  ay[7]  = stmp;
44729a97285SShri Abhyankar       stmp = ax[8];  ax[8]  = ay[8];  ay[8]  = stmp;
44829a97285SShri Abhyankar       stmp = ax[9];  ax[9]  = ay[9];  ay[9]  = stmp;
44929a97285SShri Abhyankar       stmp = ax[10]; ax[10] = ay[10]; ay[10] = stmp;
45029a97285SShri Abhyankar       stmp = ax[11]; ax[11] = ay[11]; ay[11] = stmp;
45129a97285SShri Abhyankar       stmp = ax[12]; ax[12] = ay[12]; ay[12] = stmp;
45229a97285SShri Abhyankar       stmp = ax[13]; ax[13] = ay[13]; ay[13] = stmp;
45329a97285SShri Abhyankar       stmp = ax[14]; ax[14] = ay[14]; ay[14] = stmp;
45429a97285SShri Abhyankar     }
45529a97285SShri Abhyankar   }
45629a97285SShri Abhyankar   PetscFunctionReturn(0);
45729a97285SShri Abhyankar }
458