xref: /petsc/src/mat/impls/baij/seq/dgefa2.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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        This is a combination of the Linpack routines
984643e36SBarry Smith     dgefa() and dgedi() specialized for a size of 2.
1084643e36SBarry Smith 
1184643e36SBarry Smith */
12c6db04a5SJed Brown #include <petscsys.h>
1384643e36SBarry Smith 
14*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_2(MatScalar *a, PetscReal shift, PetscBool allowzeropivot, PetscBool *zeropivotdetected) {
15690b6cddSBarry Smith   PetscInt   i__2, i__3, kp1, j, k, l, ll, i, ipvt[2], k3;
16690b6cddSBarry Smith   PetscInt   k4, j3;
17b48ee343SBarry Smith   MatScalar *aa, *ax, *ay, work[4], stmp;
18329f5518SBarry Smith   MatReal    tmp, max;
1984643e36SBarry Smith 
2084643e36SBarry Smith   PetscFunctionBegin;
21c80103daSHong Zhang   if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
22943c8ff5SShri Abhyankar   shift = .25 * shift * (1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[3]));
23c9b7c560SHong Zhang 
2484643e36SBarry Smith   /* Parameter adjustments */
2584643e36SBarry Smith   a -= 3;
2684643e36SBarry Smith 
2784643e36SBarry Smith   k   = 1;
2884643e36SBarry Smith   kp1 = k + 1;
2984643e36SBarry Smith   k3  = 2 * k;
3084643e36SBarry Smith   k4  = k3 + k;
3184643e36SBarry Smith 
32c9b7c560SHong Zhang   /* find l = pivot index */
33ed33f8a5SSatish Balay   i__2 = 3 - k;
3484643e36SBarry Smith   aa   = &a[k4];
3584643e36SBarry Smith   max  = PetscAbsScalar(aa[0]);
3684643e36SBarry Smith   l    = 1;
3784643e36SBarry Smith   for (ll = 1; ll < i__2; ll++) {
3884643e36SBarry Smith     tmp = PetscAbsScalar(aa[ll]);
39*9371c9d4SSatish Balay     if (tmp > max) {
40*9371c9d4SSatish Balay       max = tmp;
41*9371c9d4SSatish Balay       l   = ll + 1;
42*9371c9d4SSatish Balay     }
4384643e36SBarry Smith   }
4484643e36SBarry Smith   l += k - 1;
45b48ee343SBarry Smith   ipvt[k - 1] = l;
4684643e36SBarry Smith 
47943c8ff5SShri Abhyankar   if (a[l + k3] == 0.0) {
48ec1892c8SHong Zhang     if (shift == 0.0) {
49ec1892c8SHong Zhang       if (allowzeropivot) {
509566063dSJacob Faibussowitsch         PetscCall(PetscInfo(NULL, "Zero pivot, row %" PetscInt_FMT "\n", k - 1));
5170d8d27cSBarry Smith         if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
5298921bdaSJacob Faibussowitsch       } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row %" PetscInt_FMT, k - 1);
53ec1892c8SHong Zhang     } else {
54943c8ff5SShri Abhyankar       a[l + k3] = shift;
55943c8ff5SShri Abhyankar     }
56943c8ff5SShri Abhyankar   }
5784643e36SBarry Smith 
5884643e36SBarry Smith   /* interchange if necessary */
5984643e36SBarry Smith   if (l != k) {
6084643e36SBarry Smith     stmp      = a[l + k3];
6184643e36SBarry Smith     a[l + k3] = a[k4];
6284643e36SBarry Smith     a[k4]     = stmp;
6384643e36SBarry Smith   }
6484643e36SBarry Smith 
6584643e36SBarry Smith   /* compute multipliers */
6684643e36SBarry Smith   stmp = -1. / a[k4];
6784643e36SBarry Smith   i__2 = 2 - k;
6884643e36SBarry Smith   aa   = &a[1 + k4];
6926fbe8dcSKarl Rupp   for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp;
7084643e36SBarry Smith 
7184643e36SBarry Smith   /* row elimination with column indexing */
7284643e36SBarry Smith   ax = &a[k4 + 1];
7384643e36SBarry Smith   for (j = kp1; j <= 2; ++j) {
7484643e36SBarry Smith     j3   = 2 * j;
7584643e36SBarry Smith     stmp = a[l + j3];
7684643e36SBarry Smith     if (l != k) {
7784643e36SBarry Smith       a[l + j3] = a[k + j3];
7884643e36SBarry Smith       a[k + j3] = stmp;
7984643e36SBarry Smith     }
8084643e36SBarry Smith 
8184643e36SBarry Smith     i__3 = 2 - k;
8284643e36SBarry Smith     ay   = &a[1 + k + j3];
8326fbe8dcSKarl Rupp     for (ll = 0; ll < i__3; ll++) ay[ll] += stmp * ax[ll];
8484643e36SBarry Smith   }
8526fbe8dcSKarl Rupp 
86b48ee343SBarry Smith   ipvt[1] = 2;
872e92ee13SHong Zhang   if (a[6] == 0.0) {
88c0aa6a63SJacob Faibussowitsch     if (PetscLikely(allowzeropivot)) {
899566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Zero pivot, row 1\n"));
9070d8d27cSBarry Smith       if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
91c0aa6a63SJacob Faibussowitsch     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row 1");
922e92ee13SHong Zhang   }
9384643e36SBarry Smith 
94c9b7c560SHong Zhang   /* Now form the inverse */
9584643e36SBarry Smith   /* compute inverse(u) */
9684643e36SBarry Smith   for (k = 1; k <= 2; ++k) {
9784643e36SBarry Smith     k3    = 2 * k;
9884643e36SBarry Smith     k4    = k3 + k;
9984643e36SBarry Smith     a[k4] = 1.0 / a[k4];
10084643e36SBarry Smith     stmp  = -a[k4];
10184643e36SBarry Smith     i__2  = k - 1;
10284643e36SBarry Smith     aa    = &a[k3 + 1];
10384643e36SBarry Smith     for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp;
10484643e36SBarry Smith     kp1 = k + 1;
10584643e36SBarry Smith     if (2 < kp1) continue;
10684643e36SBarry Smith     ax = aa;
10784643e36SBarry Smith     for (j = kp1; j <= 2; ++j) {
10884643e36SBarry Smith       j3        = 2 * j;
10984643e36SBarry Smith       stmp      = a[k + j3];
11084643e36SBarry Smith       a[k + j3] = 0.0;
11184643e36SBarry Smith       ay        = &a[j3 + 1];
11226fbe8dcSKarl Rupp       for (ll = 0; ll < k; ll++) ay[ll] += stmp * ax[ll];
11384643e36SBarry Smith     }
11484643e36SBarry Smith   }
11584643e36SBarry Smith 
11684643e36SBarry Smith   /* form inverse(u)*inverse(l) */
11784643e36SBarry Smith   k   = 1;
11884643e36SBarry Smith   k3  = 2 * k;
11984643e36SBarry Smith   kp1 = k + 1;
12084643e36SBarry Smith   aa  = a + k3;
12184643e36SBarry Smith   for (i = kp1; i <= 2; ++i) {
122b48ee343SBarry Smith     work[i - 1] = aa[i];
12384643e36SBarry Smith     aa[i]       = 0.0;
12484643e36SBarry Smith   }
12584643e36SBarry Smith   for (j = kp1; j <= 2; ++j) {
126b48ee343SBarry Smith     stmp = work[j - 1];
12784643e36SBarry Smith     ax   = &a[2 * j + 1];
12884643e36SBarry Smith     ay   = &a[k3 + 1];
12984643e36SBarry Smith     ay[0] += stmp * ax[0];
13084643e36SBarry Smith     ay[1] += stmp * ax[1];
13184643e36SBarry Smith   }
132b48ee343SBarry Smith   l = ipvt[k - 1];
13384643e36SBarry Smith   if (l != k) {
13484643e36SBarry Smith     ax    = &a[k3 + 1];
13584643e36SBarry Smith     ay    = &a[2 * l + 1];
136*9371c9d4SSatish Balay     stmp  = ax[0];
137*9371c9d4SSatish Balay     ax[0] = ay[0];
138*9371c9d4SSatish Balay     ay[0] = stmp;
139*9371c9d4SSatish Balay     stmp  = ax[1];
140*9371c9d4SSatish Balay     ax[1] = ay[1];
141*9371c9d4SSatish Balay     ay[1] = stmp;
14284643e36SBarry Smith   }
14384643e36SBarry Smith   PetscFunctionReturn(0);
14484643e36SBarry Smith }
14584643e36SBarry Smith 
146ec1892c8SHong Zhang /* gaussian elimination with partial pivoting */
147*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_9(MatScalar *a, PetscReal shift, PetscBool allowzeropivot, PetscBool *zeropivotdetected) {
1483dfa136dSBarry Smith   PetscInt   i__2, i__3, kp1, j, k, l, ll, i, ipvt[9], kb, k3;
1493dfa136dSBarry Smith   PetscInt   k4, j3;
1503dfa136dSBarry Smith   MatScalar *aa, *ax, *ay, work[81], stmp;
1513dfa136dSBarry Smith   MatReal    tmp, max;
1523dfa136dSBarry Smith 
1533dfa136dSBarry Smith   PetscFunctionBegin;
154c80103daSHong Zhang   if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
155c80103daSHong Zhang 
1563dfa136dSBarry Smith   /* Parameter adjustments */
1573dfa136dSBarry Smith   a -= 10;
1583dfa136dSBarry Smith 
1593dfa136dSBarry Smith   for (k = 1; k <= 8; ++k) {
1603dfa136dSBarry Smith     kp1 = k + 1;
1613dfa136dSBarry Smith     k3  = 9 * k;
1623dfa136dSBarry Smith     k4  = k3 + k;
1633dfa136dSBarry Smith 
164ec1892c8SHong Zhang     /* find l = pivot index */
165c38ccd74SBarry Smith     i__2 = 10 - k;
1663dfa136dSBarry Smith     aa   = &a[k4];
1673dfa136dSBarry Smith     max  = PetscAbsScalar(aa[0]);
1683dfa136dSBarry Smith     l    = 1;
1693dfa136dSBarry Smith     for (ll = 1; ll < i__2; ll++) {
1703dfa136dSBarry Smith       tmp = PetscAbsScalar(aa[ll]);
171*9371c9d4SSatish Balay       if (tmp > max) {
172*9371c9d4SSatish Balay         max = tmp;
173*9371c9d4SSatish Balay         l   = ll + 1;
174*9371c9d4SSatish Balay       }
1753dfa136dSBarry Smith     }
1763dfa136dSBarry Smith     l += k - 1;
1773dfa136dSBarry Smith     ipvt[k - 1] = l;
1783dfa136dSBarry Smith 
179ec1892c8SHong Zhang     if (a[l + k3] == 0.0) {
180ec1892c8SHong Zhang       if (shift == 0.0) {
181ec1892c8SHong Zhang         if (allowzeropivot) {
1829566063dSJacob Faibussowitsch           PetscCall(PetscInfo(NULL, "Zero pivot, row %" PetscInt_FMT "\n", k - 1));
18370d8d27cSBarry Smith           if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
18498921bdaSJacob Faibussowitsch         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row %" PetscInt_FMT, 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) {
220c0aa6a63SJacob Faibussowitsch     if (PetscLikely(allowzeropivot)) {
2219566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Zero pivot, row 8\n"));
22270d8d27cSBarry Smith       if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
223c0aa6a63SJacob Faibussowitsch     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row 8");
2242e92ee13SHong Zhang   }
2253dfa136dSBarry Smith 
226ec1892c8SHong Zhang   /* Now form the inverse */
2273dfa136dSBarry Smith   /* compute inverse(u) */
2283dfa136dSBarry Smith   for (k = 1; k <= 9; ++k) {
2293dfa136dSBarry Smith     k3    = 9 * k;
2303dfa136dSBarry Smith     k4    = k3 + k;
2313dfa136dSBarry Smith     a[k4] = 1.0 / a[k4];
2323dfa136dSBarry Smith     stmp  = -a[k4];
2333dfa136dSBarry Smith     i__2  = k - 1;
2343dfa136dSBarry Smith     aa    = &a[k3 + 1];
2353dfa136dSBarry Smith     for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp;
2363dfa136dSBarry Smith     kp1 = k + 1;
2373dfa136dSBarry Smith     if (9 < kp1) continue;
2383dfa136dSBarry Smith     ax = aa;
2393dfa136dSBarry Smith     for (j = kp1; j <= 9; ++j) {
2403dfa136dSBarry Smith       j3        = 9 * j;
2413dfa136dSBarry Smith       stmp      = a[k + j3];
2423dfa136dSBarry Smith       a[k + j3] = 0.0;
2433dfa136dSBarry Smith       ay        = &a[j3 + 1];
24426fbe8dcSKarl Rupp       for (ll = 0; ll < k; ll++) ay[ll] += stmp * ax[ll];
2453dfa136dSBarry Smith     }
2463dfa136dSBarry Smith   }
2473dfa136dSBarry Smith 
2483dfa136dSBarry Smith   /* form inverse(u)*inverse(l) */
2493dfa136dSBarry Smith   for (kb = 1; kb <= 8; ++kb) {
2503dfa136dSBarry Smith     k   = 9 - kb;
2513dfa136dSBarry Smith     k3  = 9 * k;
2523dfa136dSBarry Smith     kp1 = k + 1;
2533dfa136dSBarry Smith     aa  = a + k3;
2543dfa136dSBarry Smith     for (i = kp1; i <= 9; ++i) {
2553dfa136dSBarry Smith       work[i - 1] = aa[i];
2563dfa136dSBarry Smith       aa[i]       = 0.0;
2573dfa136dSBarry Smith     }
2583dfa136dSBarry Smith     for (j = kp1; j <= 9; ++j) {
2593dfa136dSBarry Smith       stmp = work[j - 1];
2603dfa136dSBarry Smith       ax   = &a[9 * j + 1];
2613dfa136dSBarry Smith       ay   = &a[k3 + 1];
2623dfa136dSBarry Smith       ay[0] += stmp * ax[0];
2633dfa136dSBarry Smith       ay[1] += stmp * ax[1];
2643dfa136dSBarry Smith       ay[2] += stmp * ax[2];
2653dfa136dSBarry Smith       ay[3] += stmp * ax[3];
2663dfa136dSBarry Smith       ay[4] += stmp * ax[4];
2673dfa136dSBarry Smith       ay[5] += stmp * ax[5];
2683dfa136dSBarry Smith       ay[6] += stmp * ax[6];
2693dfa136dSBarry Smith       ay[7] += stmp * ax[7];
2703dfa136dSBarry Smith       ay[8] += stmp * ax[8];
2713dfa136dSBarry Smith     }
2723dfa136dSBarry Smith     l = ipvt[k - 1];
2733dfa136dSBarry Smith     if (l != k) {
2743dfa136dSBarry Smith       ax    = &a[k3 + 1];
2753dfa136dSBarry Smith       ay    = &a[9 * l + 1];
276*9371c9d4SSatish Balay       stmp  = ax[0];
277*9371c9d4SSatish Balay       ax[0] = ay[0];
278*9371c9d4SSatish Balay       ay[0] = stmp;
279*9371c9d4SSatish Balay       stmp  = ax[1];
280*9371c9d4SSatish Balay       ax[1] = ay[1];
281*9371c9d4SSatish Balay       ay[1] = stmp;
282*9371c9d4SSatish Balay       stmp  = ax[2];
283*9371c9d4SSatish Balay       ax[2] = ay[2];
284*9371c9d4SSatish Balay       ay[2] = stmp;
285*9371c9d4SSatish Balay       stmp  = ax[3];
286*9371c9d4SSatish Balay       ax[3] = ay[3];
287*9371c9d4SSatish Balay       ay[3] = stmp;
288*9371c9d4SSatish Balay       stmp  = ax[4];
289*9371c9d4SSatish Balay       ax[4] = ay[4];
290*9371c9d4SSatish Balay       ay[4] = stmp;
291*9371c9d4SSatish Balay       stmp  = ax[5];
292*9371c9d4SSatish Balay       ax[5] = ay[5];
293*9371c9d4SSatish Balay       ay[5] = stmp;
294*9371c9d4SSatish Balay       stmp  = ax[6];
295*9371c9d4SSatish Balay       ax[6] = ay[6];
296*9371c9d4SSatish Balay       ay[6] = stmp;
297*9371c9d4SSatish Balay       stmp  = ax[7];
298*9371c9d4SSatish Balay       ax[7] = ay[7];
299*9371c9d4SSatish Balay       ay[7] = stmp;
300*9371c9d4SSatish Balay       stmp  = ax[8];
301*9371c9d4SSatish Balay       ax[8] = ay[8];
302*9371c9d4SSatish Balay       ay[8] = stmp;
3033dfa136dSBarry Smith     }
3043dfa136dSBarry Smith   }
3053dfa136dSBarry Smith   PetscFunctionReturn(0);
3063dfa136dSBarry Smith }
30784643e36SBarry Smith 
30829a97285SShri Abhyankar /*
309ec1892c8SHong Zhang       Inverts 15 by 15 matrix using gaussian elimination with partial pivoting.
31029a97285SShri Abhyankar 
31129a97285SShri Abhyankar        Used by the sparse factorization routines in
31229a97285SShri Abhyankar      src/mat/impls/baij/seq
31329a97285SShri Abhyankar 
31429a97285SShri Abhyankar        This is a combination of the Linpack routines
31529a97285SShri Abhyankar     dgefa() and dgedi() specialized for a size of 15.
31629a97285SShri Abhyankar 
31729a97285SShri Abhyankar */
31829a97285SShri Abhyankar 
319*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_15(MatScalar *a, PetscInt *ipvt, MatScalar *work, PetscReal shift, PetscBool allowzeropivot, PetscBool *zeropivotdetected) {
320766f9fbaSBarry Smith   PetscInt   i__2, i__3, kp1, j, k, l, ll, i, kb, k3;
32129a97285SShri Abhyankar   PetscInt   k4, j3;
322766f9fbaSBarry Smith   MatScalar *aa, *ax, *ay, stmp;
32329a97285SShri Abhyankar   MatReal    tmp, max;
32429a97285SShri Abhyankar 
32529a97285SShri Abhyankar   PetscFunctionBegin;
326c80103daSHong Zhang   if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
327c80103daSHong Zhang 
32829a97285SShri Abhyankar   /* Parameter adjustments */
32929a97285SShri Abhyankar   a -= 16;
33029a97285SShri Abhyankar 
33129a97285SShri Abhyankar   for (k = 1; k <= 14; ++k) {
33229a97285SShri Abhyankar     kp1 = k + 1;
33329a97285SShri Abhyankar     k3  = 15 * k;
33429a97285SShri Abhyankar     k4  = k3 + k;
33529a97285SShri Abhyankar 
336ec1892c8SHong Zhang     /* find l = pivot index */
33729a97285SShri Abhyankar     i__2 = 16 - k;
33829a97285SShri Abhyankar     aa   = &a[k4];
33929a97285SShri Abhyankar     max  = PetscAbsScalar(aa[0]);
34029a97285SShri Abhyankar     l    = 1;
34129a97285SShri Abhyankar     for (ll = 1; ll < i__2; ll++) {
34229a97285SShri Abhyankar       tmp = PetscAbsScalar(aa[ll]);
343*9371c9d4SSatish Balay       if (tmp > max) {
344*9371c9d4SSatish Balay         max = tmp;
345*9371c9d4SSatish Balay         l   = ll + 1;
346*9371c9d4SSatish Balay       }
34729a97285SShri Abhyankar     }
34829a97285SShri Abhyankar     l += k - 1;
34929a97285SShri Abhyankar     ipvt[k - 1] = l;
35029a97285SShri Abhyankar 
351ec1892c8SHong Zhang     if (a[l + k3] == 0.0) {
352ec1892c8SHong Zhang       if (shift == 0.0) {
353ec1892c8SHong Zhang         if (allowzeropivot) {
3549566063dSJacob Faibussowitsch           PetscCall(PetscInfo(NULL, "Zero pivot, row %" PetscInt_FMT "\n", k - 1));
35570d8d27cSBarry Smith           if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
35698921bdaSJacob Faibussowitsch         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row %" PetscInt_FMT, k - 1);
357ec1892c8SHong Zhang       } else {
358ec1892c8SHong Zhang         a[l + k3] = shift;
359ec1892c8SHong Zhang       }
360ec1892c8SHong Zhang     }
36129a97285SShri Abhyankar 
36229a97285SShri Abhyankar     /* interchange if necessary */
36329a97285SShri Abhyankar     if (l != k) {
36429a97285SShri Abhyankar       stmp      = a[l + k3];
36529a97285SShri Abhyankar       a[l + k3] = a[k4];
36629a97285SShri Abhyankar       a[k4]     = stmp;
36729a97285SShri Abhyankar     }
36829a97285SShri Abhyankar 
36929a97285SShri Abhyankar     /* compute multipliers */
37029a97285SShri Abhyankar     stmp = -1. / a[k4];
37129a97285SShri Abhyankar     i__2 = 15 - k;
37229a97285SShri Abhyankar     aa   = &a[1 + k4];
37326fbe8dcSKarl Rupp     for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp;
37429a97285SShri Abhyankar 
37529a97285SShri Abhyankar     /* row elimination with column indexing */
37629a97285SShri Abhyankar     ax = &a[k4 + 1];
37729a97285SShri Abhyankar     for (j = kp1; j <= 15; ++j) {
37829a97285SShri Abhyankar       j3   = 15 * j;
37929a97285SShri Abhyankar       stmp = a[l + j3];
38029a97285SShri Abhyankar       if (l != k) {
38129a97285SShri Abhyankar         a[l + j3] = a[k + j3];
38229a97285SShri Abhyankar         a[k + j3] = stmp;
38329a97285SShri Abhyankar       }
38429a97285SShri Abhyankar 
38529a97285SShri Abhyankar       i__3 = 15 - k;
38629a97285SShri Abhyankar       ay   = &a[1 + k + j3];
38726fbe8dcSKarl Rupp       for (ll = 0; ll < i__3; ll++) ay[ll] += stmp * ax[ll];
38829a97285SShri Abhyankar     }
38929a97285SShri Abhyankar   }
39029a97285SShri Abhyankar   ipvt[14] = 15;
391922032d7SHong Zhang   if (a[240] == 0.0) {
392c0aa6a63SJacob Faibussowitsch     if (PetscLikely(allowzeropivot)) {
3939566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Zero pivot, row 14\n"));
39470d8d27cSBarry Smith       if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
395c0aa6a63SJacob Faibussowitsch     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row 14");
396922032d7SHong Zhang   }
39729a97285SShri Abhyankar 
398ec1892c8SHong Zhang   /* Now form the inverse */
39929a97285SShri Abhyankar   /* compute inverse(u) */
40029a97285SShri Abhyankar   for (k = 1; k <= 15; ++k) {
40129a97285SShri Abhyankar     k3    = 15 * k;
40229a97285SShri Abhyankar     k4    = k3 + k;
40329a97285SShri Abhyankar     a[k4] = 1.0 / a[k4];
40429a97285SShri Abhyankar     stmp  = -a[k4];
40529a97285SShri Abhyankar     i__2  = k - 1;
40629a97285SShri Abhyankar     aa    = &a[k3 + 1];
40729a97285SShri Abhyankar     for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp;
40829a97285SShri Abhyankar     kp1 = k + 1;
40929a97285SShri Abhyankar     if (15 < kp1) continue;
41029a97285SShri Abhyankar     ax = aa;
41129a97285SShri Abhyankar     for (j = kp1; j <= 15; ++j) {
41229a97285SShri Abhyankar       j3        = 15 * j;
41329a97285SShri Abhyankar       stmp      = a[k + j3];
41429a97285SShri Abhyankar       a[k + j3] = 0.0;
41529a97285SShri Abhyankar       ay        = &a[j3 + 1];
41626fbe8dcSKarl Rupp       for (ll = 0; ll < k; ll++) ay[ll] += stmp * ax[ll];
41729a97285SShri Abhyankar     }
41829a97285SShri Abhyankar   }
41929a97285SShri Abhyankar 
42029a97285SShri Abhyankar   /* form inverse(u)*inverse(l) */
42129a97285SShri Abhyankar   for (kb = 1; kb <= 14; ++kb) {
42229a97285SShri Abhyankar     k   = 15 - kb;
42329a97285SShri Abhyankar     k3  = 15 * k;
42429a97285SShri Abhyankar     kp1 = k + 1;
42529a97285SShri Abhyankar     aa  = a + k3;
42629a97285SShri Abhyankar     for (i = kp1; i <= 15; ++i) {
42729a97285SShri Abhyankar       work[i - 1] = aa[i];
42829a97285SShri Abhyankar       aa[i]       = 0.0;
42929a97285SShri Abhyankar     }
43029a97285SShri Abhyankar     for (j = kp1; j <= 15; ++j) {
43129a97285SShri Abhyankar       stmp = work[j - 1];
43229a97285SShri Abhyankar       ax   = &a[15 * j + 1];
43329a97285SShri Abhyankar       ay   = &a[k3 + 1];
43429a97285SShri Abhyankar       ay[0] += stmp * ax[0];
43529a97285SShri Abhyankar       ay[1] += stmp * ax[1];
43629a97285SShri Abhyankar       ay[2] += stmp * ax[2];
43729a97285SShri Abhyankar       ay[3] += stmp * ax[3];
43829a97285SShri Abhyankar       ay[4] += stmp * ax[4];
43929a97285SShri Abhyankar       ay[5] += stmp * ax[5];
44029a97285SShri Abhyankar       ay[6] += stmp * ax[6];
44129a97285SShri Abhyankar       ay[7] += stmp * ax[7];
44229a97285SShri Abhyankar       ay[8] += stmp * ax[8];
44329a97285SShri Abhyankar       ay[9] += stmp * ax[9];
44429a97285SShri Abhyankar       ay[10] += stmp * ax[10];
44529a97285SShri Abhyankar       ay[11] += stmp * ax[11];
44629a97285SShri Abhyankar       ay[12] += stmp * ax[12];
44729a97285SShri Abhyankar       ay[13] += stmp * ax[13];
44829a97285SShri Abhyankar       ay[14] += stmp * ax[14];
44929a97285SShri Abhyankar     }
45029a97285SShri Abhyankar     l = ipvt[k - 1];
45129a97285SShri Abhyankar     if (l != k) {
45229a97285SShri Abhyankar       ax     = &a[k3 + 1];
45329a97285SShri Abhyankar       ay     = &a[15 * l + 1];
454*9371c9d4SSatish Balay       stmp   = ax[0];
455*9371c9d4SSatish Balay       ax[0]  = ay[0];
456*9371c9d4SSatish Balay       ay[0]  = stmp;
457*9371c9d4SSatish Balay       stmp   = ax[1];
458*9371c9d4SSatish Balay       ax[1]  = ay[1];
459*9371c9d4SSatish Balay       ay[1]  = stmp;
460*9371c9d4SSatish Balay       stmp   = ax[2];
461*9371c9d4SSatish Balay       ax[2]  = ay[2];
462*9371c9d4SSatish Balay       ay[2]  = stmp;
463*9371c9d4SSatish Balay       stmp   = ax[3];
464*9371c9d4SSatish Balay       ax[3]  = ay[3];
465*9371c9d4SSatish Balay       ay[3]  = stmp;
466*9371c9d4SSatish Balay       stmp   = ax[4];
467*9371c9d4SSatish Balay       ax[4]  = ay[4];
468*9371c9d4SSatish Balay       ay[4]  = stmp;
469*9371c9d4SSatish Balay       stmp   = ax[5];
470*9371c9d4SSatish Balay       ax[5]  = ay[5];
471*9371c9d4SSatish Balay       ay[5]  = stmp;
472*9371c9d4SSatish Balay       stmp   = ax[6];
473*9371c9d4SSatish Balay       ax[6]  = ay[6];
474*9371c9d4SSatish Balay       ay[6]  = stmp;
475*9371c9d4SSatish Balay       stmp   = ax[7];
476*9371c9d4SSatish Balay       ax[7]  = ay[7];
477*9371c9d4SSatish Balay       ay[7]  = stmp;
478*9371c9d4SSatish Balay       stmp   = ax[8];
479*9371c9d4SSatish Balay       ax[8]  = ay[8];
480*9371c9d4SSatish Balay       ay[8]  = stmp;
481*9371c9d4SSatish Balay       stmp   = ax[9];
482*9371c9d4SSatish Balay       ax[9]  = ay[9];
483*9371c9d4SSatish Balay       ay[9]  = stmp;
484*9371c9d4SSatish Balay       stmp   = ax[10];
485*9371c9d4SSatish Balay       ax[10] = ay[10];
486*9371c9d4SSatish Balay       ay[10] = stmp;
487*9371c9d4SSatish Balay       stmp   = ax[11];
488*9371c9d4SSatish Balay       ax[11] = ay[11];
489*9371c9d4SSatish Balay       ay[11] = stmp;
490*9371c9d4SSatish Balay       stmp   = ax[12];
491*9371c9d4SSatish Balay       ax[12] = ay[12];
492*9371c9d4SSatish Balay       ay[12] = stmp;
493*9371c9d4SSatish Balay       stmp   = ax[13];
494*9371c9d4SSatish Balay       ax[13] = ay[13];
495*9371c9d4SSatish Balay       ay[13] = stmp;
496*9371c9d4SSatish Balay       stmp   = ax[14];
497*9371c9d4SSatish Balay       ax[14] = ay[14];
498*9371c9d4SSatish Balay       ay[14] = stmp;
49929a97285SShri Abhyankar     }
50029a97285SShri Abhyankar   }
50129a97285SShri Abhyankar   PetscFunctionReturn(0);
50229a97285SShri Abhyankar }
503