184643e36SBarry Smith /*
2c9b7c560SHong Zhang Inverts 2 by 2 matrix using gaussian elimination with partial pivoting.
384643e36SBarry Smith
484643e36SBarry Smith Used by the sparse factorization routines in
5dd882469SBarry Smith src/mat/impls/baij/seq
684643e36SBarry Smith
784643e36SBarry Smith This is a combination of the Linpack routines
884643e36SBarry Smith dgefa() and dgedi() specialized for a size of 2.
984643e36SBarry Smith
1084643e36SBarry Smith */
11c6db04a5SJed Brown #include <petscsys.h>
12d6acfc2dSPierre Jolivet #include <petsc/private/kernels/blockinvert.h>
1384643e36SBarry Smith
PetscKernel_A_gets_inverse_A_2(MatScalar * a,PetscReal shift,PetscBool allowzeropivot,PetscBool * zeropivotdetected)14d6acfc2dSPierre Jolivet PetscErrorCode PetscKernel_A_gets_inverse_A_2(MatScalar *a, PetscReal shift, PetscBool allowzeropivot, PetscBool *zeropivotdetected)
15d71ae5a4SJacob Faibussowitsch {
16690b6cddSBarry Smith PetscInt i__2, i__3, kp1, j, k, l, ll, i, ipvt[2], k3;
17690b6cddSBarry Smith PetscInt k4, j3;
18b48ee343SBarry Smith MatScalar *aa, *ax, *ay, work[4], stmp;
19329f5518SBarry Smith MatReal tmp, max;
2084643e36SBarry Smith
2184643e36SBarry Smith PetscFunctionBegin;
22c80103daSHong Zhang if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
23943c8ff5SShri Abhyankar shift = .25 * shift * (1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[3]));
24c9b7c560SHong Zhang
2584643e36SBarry Smith /* Parameter adjustments */
2684643e36SBarry Smith a -= 3;
2784643e36SBarry Smith
2884643e36SBarry Smith k = 1;
2984643e36SBarry Smith kp1 = k + 1;
3084643e36SBarry Smith k3 = 2 * k;
3184643e36SBarry Smith k4 = k3 + k;
3284643e36SBarry Smith
33c9b7c560SHong Zhang /* find l = pivot index */
34ed33f8a5SSatish Balay i__2 = 3 - k;
3584643e36SBarry Smith aa = &a[k4];
3684643e36SBarry Smith max = PetscAbsScalar(aa[0]);
3784643e36SBarry Smith l = 1;
3884643e36SBarry Smith for (ll = 1; ll < i__2; ll++) {
3984643e36SBarry Smith tmp = PetscAbsScalar(aa[ll]);
409371c9d4SSatish Balay if (tmp > max) {
419371c9d4SSatish Balay max = tmp;
429371c9d4SSatish Balay l = ll + 1;
439371c9d4SSatish Balay }
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) {
50*966bd95aSPierre Jolivet PetscCheck(allowzeropivot, PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row %" PetscInt_FMT, k - 1);
519566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Zero pivot, row %" PetscInt_FMT "\n", k - 1));
5270d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
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) {
88*966bd95aSPierre Jolivet PetscCheck(PetscLikely(allowzeropivot), PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row 1");
899566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Zero pivot, row 1\n"));
9070d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
912e92ee13SHong Zhang }
9284643e36SBarry Smith
93c9b7c560SHong Zhang /* Now form the inverse */
9484643e36SBarry Smith /* compute inverse(u) */
9584643e36SBarry Smith for (k = 1; k <= 2; ++k) {
9684643e36SBarry Smith k3 = 2 * k;
9784643e36SBarry Smith k4 = k3 + k;
9884643e36SBarry Smith a[k4] = 1.0 / a[k4];
9984643e36SBarry Smith stmp = -a[k4];
10084643e36SBarry Smith i__2 = k - 1;
10184643e36SBarry Smith aa = &a[k3 + 1];
10284643e36SBarry Smith for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp;
10384643e36SBarry Smith kp1 = k + 1;
10484643e36SBarry Smith if (2 < kp1) continue;
10584643e36SBarry Smith ax = aa;
10684643e36SBarry Smith for (j = kp1; j <= 2; ++j) {
10784643e36SBarry Smith j3 = 2 * j;
10884643e36SBarry Smith stmp = a[k + j3];
10984643e36SBarry Smith a[k + j3] = 0.0;
11084643e36SBarry Smith ay = &a[j3 + 1];
11126fbe8dcSKarl Rupp for (ll = 0; ll < k; ll++) ay[ll] += stmp * ax[ll];
11284643e36SBarry Smith }
11384643e36SBarry Smith }
11484643e36SBarry Smith
11584643e36SBarry Smith /* form inverse(u)*inverse(l) */
11684643e36SBarry Smith k = 1;
11784643e36SBarry Smith k3 = 2 * k;
11884643e36SBarry Smith kp1 = k + 1;
11984643e36SBarry Smith aa = a + k3;
12084643e36SBarry Smith for (i = kp1; i <= 2; ++i) {
121b48ee343SBarry Smith work[i - 1] = aa[i];
12284643e36SBarry Smith aa[i] = 0.0;
12384643e36SBarry Smith }
12484643e36SBarry Smith for (j = kp1; j <= 2; ++j) {
125b48ee343SBarry Smith stmp = work[j - 1];
12684643e36SBarry Smith ax = &a[2 * j + 1];
12784643e36SBarry Smith ay = &a[k3 + 1];
12884643e36SBarry Smith ay[0] += stmp * ax[0];
12984643e36SBarry Smith ay[1] += stmp * ax[1];
13084643e36SBarry Smith }
131b48ee343SBarry Smith l = ipvt[k - 1];
13284643e36SBarry Smith if (l != k) {
13384643e36SBarry Smith ax = &a[k3 + 1];
13484643e36SBarry Smith ay = &a[2 * l + 1];
1359371c9d4SSatish Balay stmp = ax[0];
1369371c9d4SSatish Balay ax[0] = ay[0];
1379371c9d4SSatish Balay ay[0] = stmp;
1389371c9d4SSatish Balay stmp = ax[1];
1399371c9d4SSatish Balay ax[1] = ay[1];
1409371c9d4SSatish Balay ay[1] = stmp;
14184643e36SBarry Smith }
1423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
14384643e36SBarry Smith }
14484643e36SBarry Smith
145d6acfc2dSPierre Jolivet /* Gaussian elimination with partial pivoting */
PetscKernel_A_gets_inverse_A_9(MatScalar * a,PetscReal shift,PetscBool allowzeropivot,PetscBool * zeropivotdetected)146d6acfc2dSPierre Jolivet PetscErrorCode PetscKernel_A_gets_inverse_A_9(MatScalar *a, PetscReal shift, PetscBool allowzeropivot, PetscBool *zeropivotdetected)
147d71ae5a4SJacob Faibussowitsch {
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]);
1719371c9d4SSatish Balay if (tmp > max) {
1729371c9d4SSatish Balay max = tmp;
1739371c9d4SSatish Balay l = ll + 1;
1749371c9d4SSatish 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) {
181*966bd95aSPierre Jolivet PetscCheck(allowzeropivot, PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row %" PetscInt_FMT, k - 1);
1829566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Zero pivot, row %" PetscInt_FMT "\n", k - 1));
18370d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
184ec1892c8SHong Zhang } else {
185ec1892c8SHong Zhang a[l + k3] = shift;
186ec1892c8SHong Zhang }
187ec1892c8SHong Zhang }
1883dfa136dSBarry Smith
1893dfa136dSBarry Smith /* interchange if necessary */
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 stmp = -1. / a[k4];
1983dfa136dSBarry Smith i__2 = 9 - k;
1993dfa136dSBarry Smith aa = &a[1 + k4];
20026fbe8dcSKarl Rupp for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp;
2013dfa136dSBarry Smith
2023dfa136dSBarry Smith /* row elimination with column indexing */
2033dfa136dSBarry Smith ax = &a[k4 + 1];
2043dfa136dSBarry Smith for (j = kp1; j <= 9; ++j) {
2053dfa136dSBarry Smith j3 = 9 * j;
2063dfa136dSBarry Smith stmp = a[l + j3];
2073dfa136dSBarry Smith if (l != k) {
2083dfa136dSBarry Smith a[l + j3] = a[k + j3];
2093dfa136dSBarry Smith a[k + j3] = stmp;
2103dfa136dSBarry Smith }
2113dfa136dSBarry Smith
2123dfa136dSBarry Smith i__3 = 9 - k;
2133dfa136dSBarry Smith ay = &a[1 + k + j3];
21426fbe8dcSKarl Rupp for (ll = 0; ll < i__3; ll++) ay[ll] += stmp * ax[ll];
2153dfa136dSBarry Smith }
2163dfa136dSBarry Smith }
2173dfa136dSBarry Smith ipvt[8] = 9;
2182e92ee13SHong Zhang if (a[90] == 0.0) {
219*966bd95aSPierre Jolivet PetscCheck(PetscLikely(allowzeropivot), PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row 8");
2209566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Zero pivot, row 8\n"));
22170d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
2222e92ee13SHong Zhang }
2233dfa136dSBarry Smith
224ec1892c8SHong Zhang /* Now form the inverse */
2253dfa136dSBarry Smith /* compute inverse(u) */
2263dfa136dSBarry Smith for (k = 1; k <= 9; ++k) {
2273dfa136dSBarry Smith k3 = 9 * k;
2283dfa136dSBarry Smith k4 = k3 + k;
2293dfa136dSBarry Smith a[k4] = 1.0 / a[k4];
2303dfa136dSBarry Smith stmp = -a[k4];
2313dfa136dSBarry Smith i__2 = k - 1;
2323dfa136dSBarry Smith aa = &a[k3 + 1];
2333dfa136dSBarry Smith for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp;
2343dfa136dSBarry Smith kp1 = k + 1;
2353dfa136dSBarry Smith if (9 < kp1) continue;
2363dfa136dSBarry Smith ax = aa;
2373dfa136dSBarry Smith for (j = kp1; j <= 9; ++j) {
2383dfa136dSBarry Smith j3 = 9 * j;
2393dfa136dSBarry Smith stmp = a[k + j3];
2403dfa136dSBarry Smith a[k + j3] = 0.0;
2413dfa136dSBarry Smith ay = &a[j3 + 1];
24226fbe8dcSKarl Rupp for (ll = 0; ll < k; ll++) ay[ll] += stmp * ax[ll];
2433dfa136dSBarry Smith }
2443dfa136dSBarry Smith }
2453dfa136dSBarry Smith
2463dfa136dSBarry Smith /* form inverse(u)*inverse(l) */
2473dfa136dSBarry Smith for (kb = 1; kb <= 8; ++kb) {
2483dfa136dSBarry Smith k = 9 - kb;
2493dfa136dSBarry Smith k3 = 9 * k;
2503dfa136dSBarry Smith kp1 = k + 1;
2513dfa136dSBarry Smith aa = a + k3;
2523dfa136dSBarry Smith for (i = kp1; i <= 9; ++i) {
2533dfa136dSBarry Smith work[i - 1] = aa[i];
2543dfa136dSBarry Smith aa[i] = 0.0;
2553dfa136dSBarry Smith }
2563dfa136dSBarry Smith for (j = kp1; j <= 9; ++j) {
2573dfa136dSBarry Smith stmp = work[j - 1];
2583dfa136dSBarry Smith ax = &a[9 * j + 1];
2593dfa136dSBarry Smith ay = &a[k3 + 1];
2603dfa136dSBarry Smith ay[0] += stmp * ax[0];
2613dfa136dSBarry Smith ay[1] += stmp * ax[1];
2623dfa136dSBarry Smith ay[2] += stmp * ax[2];
2633dfa136dSBarry Smith ay[3] += stmp * ax[3];
2643dfa136dSBarry Smith ay[4] += stmp * ax[4];
2653dfa136dSBarry Smith ay[5] += stmp * ax[5];
2663dfa136dSBarry Smith ay[6] += stmp * ax[6];
2673dfa136dSBarry Smith ay[7] += stmp * ax[7];
2683dfa136dSBarry Smith ay[8] += stmp * ax[8];
2693dfa136dSBarry Smith }
2703dfa136dSBarry Smith l = ipvt[k - 1];
2713dfa136dSBarry Smith if (l != k) {
2723dfa136dSBarry Smith ax = &a[k3 + 1];
2733dfa136dSBarry Smith ay = &a[9 * l + 1];
2749371c9d4SSatish Balay stmp = ax[0];
2759371c9d4SSatish Balay ax[0] = ay[0];
2769371c9d4SSatish Balay ay[0] = stmp;
2779371c9d4SSatish Balay stmp = ax[1];
2789371c9d4SSatish Balay ax[1] = ay[1];
2799371c9d4SSatish Balay ay[1] = stmp;
2809371c9d4SSatish Balay stmp = ax[2];
2819371c9d4SSatish Balay ax[2] = ay[2];
2829371c9d4SSatish Balay ay[2] = stmp;
2839371c9d4SSatish Balay stmp = ax[3];
2849371c9d4SSatish Balay ax[3] = ay[3];
2859371c9d4SSatish Balay ay[3] = stmp;
2869371c9d4SSatish Balay stmp = ax[4];
2879371c9d4SSatish Balay ax[4] = ay[4];
2889371c9d4SSatish Balay ay[4] = stmp;
2899371c9d4SSatish Balay stmp = ax[5];
2909371c9d4SSatish Balay ax[5] = ay[5];
2919371c9d4SSatish Balay ay[5] = stmp;
2929371c9d4SSatish Balay stmp = ax[6];
2939371c9d4SSatish Balay ax[6] = ay[6];
2949371c9d4SSatish Balay ay[6] = stmp;
2959371c9d4SSatish Balay stmp = ax[7];
2969371c9d4SSatish Balay ax[7] = ay[7];
2979371c9d4SSatish Balay ay[7] = stmp;
2989371c9d4SSatish Balay stmp = ax[8];
2999371c9d4SSatish Balay ax[8] = ay[8];
3009371c9d4SSatish Balay ay[8] = stmp;
3013dfa136dSBarry Smith }
3023dfa136dSBarry Smith }
3033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3043dfa136dSBarry Smith }
30584643e36SBarry Smith
30629a97285SShri Abhyankar /*
307ec1892c8SHong Zhang Inverts 15 by 15 matrix using gaussian elimination with partial pivoting.
30829a97285SShri Abhyankar
30929a97285SShri Abhyankar Used by the sparse factorization routines in
31029a97285SShri Abhyankar src/mat/impls/baij/seq
31129a97285SShri Abhyankar
31229a97285SShri Abhyankar This is a combination of the Linpack routines
31329a97285SShri Abhyankar dgefa() and dgedi() specialized for a size of 15.
31429a97285SShri Abhyankar
31529a97285SShri Abhyankar */
31629a97285SShri Abhyankar
PetscKernel_A_gets_inverse_A_15(MatScalar * a,PetscInt * ipvt,MatScalar * work,PetscReal shift,PetscBool allowzeropivot,PetscBool * zeropivotdetected)317d6acfc2dSPierre Jolivet PetscErrorCode PetscKernel_A_gets_inverse_A_15(MatScalar *a, PetscInt *ipvt, MatScalar *work, PetscReal shift, PetscBool allowzeropivot, PetscBool *zeropivotdetected)
318d71ae5a4SJacob Faibussowitsch {
319766f9fbaSBarry Smith PetscInt i__2, i__3, kp1, j, k, l, ll, i, kb, k3;
32029a97285SShri Abhyankar PetscInt k4, j3;
321766f9fbaSBarry Smith MatScalar *aa, *ax, *ay, stmp;
32229a97285SShri Abhyankar MatReal tmp, max;
32329a97285SShri Abhyankar
32429a97285SShri Abhyankar PetscFunctionBegin;
325c80103daSHong Zhang if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
326c80103daSHong Zhang
32729a97285SShri Abhyankar /* Parameter adjustments */
32829a97285SShri Abhyankar a -= 16;
32929a97285SShri Abhyankar
33029a97285SShri Abhyankar for (k = 1; k <= 14; ++k) {
33129a97285SShri Abhyankar kp1 = k + 1;
33229a97285SShri Abhyankar k3 = 15 * k;
33329a97285SShri Abhyankar k4 = k3 + k;
33429a97285SShri Abhyankar
335ec1892c8SHong Zhang /* find l = pivot index */
33629a97285SShri Abhyankar i__2 = 16 - k;
33729a97285SShri Abhyankar aa = &a[k4];
33829a97285SShri Abhyankar max = PetscAbsScalar(aa[0]);
33929a97285SShri Abhyankar l = 1;
34029a97285SShri Abhyankar for (ll = 1; ll < i__2; ll++) {
34129a97285SShri Abhyankar tmp = PetscAbsScalar(aa[ll]);
3429371c9d4SSatish Balay if (tmp > max) {
3439371c9d4SSatish Balay max = tmp;
3449371c9d4SSatish Balay l = ll + 1;
3459371c9d4SSatish Balay }
34629a97285SShri Abhyankar }
34729a97285SShri Abhyankar l += k - 1;
34829a97285SShri Abhyankar ipvt[k - 1] = l;
34929a97285SShri Abhyankar
350ec1892c8SHong Zhang if (a[l + k3] == 0.0) {
351ec1892c8SHong Zhang if (shift == 0.0) {
352*966bd95aSPierre Jolivet PetscCheck(allowzeropivot, PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row %" PetscInt_FMT, k - 1);
3539566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Zero pivot, row %" PetscInt_FMT "\n", k - 1));
35470d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
355ec1892c8SHong Zhang } else {
356ec1892c8SHong Zhang a[l + k3] = shift;
357ec1892c8SHong Zhang }
358ec1892c8SHong Zhang }
35929a97285SShri Abhyankar
36029a97285SShri Abhyankar /* interchange if necessary */
36129a97285SShri Abhyankar if (l != k) {
36229a97285SShri Abhyankar stmp = a[l + k3];
36329a97285SShri Abhyankar a[l + k3] = a[k4];
36429a97285SShri Abhyankar a[k4] = stmp;
36529a97285SShri Abhyankar }
36629a97285SShri Abhyankar
36729a97285SShri Abhyankar /* compute multipliers */
36829a97285SShri Abhyankar stmp = -1. / a[k4];
36929a97285SShri Abhyankar i__2 = 15 - k;
37029a97285SShri Abhyankar aa = &a[1 + k4];
37126fbe8dcSKarl Rupp for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp;
37229a97285SShri Abhyankar
37329a97285SShri Abhyankar /* row elimination with column indexing */
37429a97285SShri Abhyankar ax = &a[k4 + 1];
37529a97285SShri Abhyankar for (j = kp1; j <= 15; ++j) {
37629a97285SShri Abhyankar j3 = 15 * j;
37729a97285SShri Abhyankar stmp = a[l + j3];
37829a97285SShri Abhyankar if (l != k) {
37929a97285SShri Abhyankar a[l + j3] = a[k + j3];
38029a97285SShri Abhyankar a[k + j3] = stmp;
38129a97285SShri Abhyankar }
38229a97285SShri Abhyankar
38329a97285SShri Abhyankar i__3 = 15 - k;
38429a97285SShri Abhyankar ay = &a[1 + k + j3];
38526fbe8dcSKarl Rupp for (ll = 0; ll < i__3; ll++) ay[ll] += stmp * ax[ll];
38629a97285SShri Abhyankar }
38729a97285SShri Abhyankar }
38829a97285SShri Abhyankar ipvt[14] = 15;
389922032d7SHong Zhang if (a[240] == 0.0) {
390*966bd95aSPierre Jolivet PetscCheck(PetscLikely(allowzeropivot), PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row 14");
3919566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Zero pivot, row 14\n"));
39270d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
393922032d7SHong Zhang }
39429a97285SShri Abhyankar
395ec1892c8SHong Zhang /* Now form the inverse */
39629a97285SShri Abhyankar /* compute inverse(u) */
39729a97285SShri Abhyankar for (k = 1; k <= 15; ++k) {
39829a97285SShri Abhyankar k3 = 15 * k;
39929a97285SShri Abhyankar k4 = k3 + k;
40029a97285SShri Abhyankar a[k4] = 1.0 / a[k4];
40129a97285SShri Abhyankar stmp = -a[k4];
40229a97285SShri Abhyankar i__2 = k - 1;
40329a97285SShri Abhyankar aa = &a[k3 + 1];
40429a97285SShri Abhyankar for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp;
40529a97285SShri Abhyankar kp1 = k + 1;
40629a97285SShri Abhyankar if (15 < kp1) continue;
40729a97285SShri Abhyankar ax = aa;
40829a97285SShri Abhyankar for (j = kp1; j <= 15; ++j) {
40929a97285SShri Abhyankar j3 = 15 * j;
41029a97285SShri Abhyankar stmp = a[k + j3];
41129a97285SShri Abhyankar a[k + j3] = 0.0;
41229a97285SShri Abhyankar ay = &a[j3 + 1];
41326fbe8dcSKarl Rupp for (ll = 0; ll < k; ll++) ay[ll] += stmp * ax[ll];
41429a97285SShri Abhyankar }
41529a97285SShri Abhyankar }
41629a97285SShri Abhyankar
41729a97285SShri Abhyankar /* form inverse(u)*inverse(l) */
41829a97285SShri Abhyankar for (kb = 1; kb <= 14; ++kb) {
41929a97285SShri Abhyankar k = 15 - kb;
42029a97285SShri Abhyankar k3 = 15 * k;
42129a97285SShri Abhyankar kp1 = k + 1;
42229a97285SShri Abhyankar aa = a + k3;
42329a97285SShri Abhyankar for (i = kp1; i <= 15; ++i) {
42429a97285SShri Abhyankar work[i - 1] = aa[i];
42529a97285SShri Abhyankar aa[i] = 0.0;
42629a97285SShri Abhyankar }
42729a97285SShri Abhyankar for (j = kp1; j <= 15; ++j) {
42829a97285SShri Abhyankar stmp = work[j - 1];
42929a97285SShri Abhyankar ax = &a[15 * j + 1];
43029a97285SShri Abhyankar ay = &a[k3 + 1];
43129a97285SShri Abhyankar ay[0] += stmp * ax[0];
43229a97285SShri Abhyankar ay[1] += stmp * ax[1];
43329a97285SShri Abhyankar ay[2] += stmp * ax[2];
43429a97285SShri Abhyankar ay[3] += stmp * ax[3];
43529a97285SShri Abhyankar ay[4] += stmp * ax[4];
43629a97285SShri Abhyankar ay[5] += stmp * ax[5];
43729a97285SShri Abhyankar ay[6] += stmp * ax[6];
43829a97285SShri Abhyankar ay[7] += stmp * ax[7];
43929a97285SShri Abhyankar ay[8] += stmp * ax[8];
44029a97285SShri Abhyankar ay[9] += stmp * ax[9];
44129a97285SShri Abhyankar ay[10] += stmp * ax[10];
44229a97285SShri Abhyankar ay[11] += stmp * ax[11];
44329a97285SShri Abhyankar ay[12] += stmp * ax[12];
44429a97285SShri Abhyankar ay[13] += stmp * ax[13];
44529a97285SShri Abhyankar ay[14] += stmp * ax[14];
44629a97285SShri Abhyankar }
44729a97285SShri Abhyankar l = ipvt[k - 1];
44829a97285SShri Abhyankar if (l != k) {
44929a97285SShri Abhyankar ax = &a[k3 + 1];
45029a97285SShri Abhyankar ay = &a[15 * l + 1];
4519371c9d4SSatish Balay stmp = ax[0];
4529371c9d4SSatish Balay ax[0] = ay[0];
4539371c9d4SSatish Balay ay[0] = stmp;
4549371c9d4SSatish Balay stmp = ax[1];
4559371c9d4SSatish Balay ax[1] = ay[1];
4569371c9d4SSatish Balay ay[1] = stmp;
4579371c9d4SSatish Balay stmp = ax[2];
4589371c9d4SSatish Balay ax[2] = ay[2];
4599371c9d4SSatish Balay ay[2] = stmp;
4609371c9d4SSatish Balay stmp = ax[3];
4619371c9d4SSatish Balay ax[3] = ay[3];
4629371c9d4SSatish Balay ay[3] = stmp;
4639371c9d4SSatish Balay stmp = ax[4];
4649371c9d4SSatish Balay ax[4] = ay[4];
4659371c9d4SSatish Balay ay[4] = stmp;
4669371c9d4SSatish Balay stmp = ax[5];
4679371c9d4SSatish Balay ax[5] = ay[5];
4689371c9d4SSatish Balay ay[5] = stmp;
4699371c9d4SSatish Balay stmp = ax[6];
4709371c9d4SSatish Balay ax[6] = ay[6];
4719371c9d4SSatish Balay ay[6] = stmp;
4729371c9d4SSatish Balay stmp = ax[7];
4739371c9d4SSatish Balay ax[7] = ay[7];
4749371c9d4SSatish Balay ay[7] = stmp;
4759371c9d4SSatish Balay stmp = ax[8];
4769371c9d4SSatish Balay ax[8] = ay[8];
4779371c9d4SSatish Balay ay[8] = stmp;
4789371c9d4SSatish Balay stmp = ax[9];
4799371c9d4SSatish Balay ax[9] = ay[9];
4809371c9d4SSatish Balay ay[9] = stmp;
4819371c9d4SSatish Balay stmp = ax[10];
4829371c9d4SSatish Balay ax[10] = ay[10];
4839371c9d4SSatish Balay ay[10] = stmp;
4849371c9d4SSatish Balay stmp = ax[11];
4859371c9d4SSatish Balay ax[11] = ay[11];
4869371c9d4SSatish Balay ay[11] = stmp;
4879371c9d4SSatish Balay stmp = ax[12];
4889371c9d4SSatish Balay ax[12] = ay[12];
4899371c9d4SSatish Balay ay[12] = stmp;
4909371c9d4SSatish Balay stmp = ax[13];
4919371c9d4SSatish Balay ax[13] = ay[13];
4929371c9d4SSatish Balay ay[13] = stmp;
4939371c9d4SSatish Balay stmp = ax[14];
4949371c9d4SSatish Balay ax[14] = ay[14];
4959371c9d4SSatish Balay ay[14] = stmp;
49629a97285SShri Abhyankar }
49729a97285SShri Abhyankar }
4983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
49929a97285SShri Abhyankar }
500