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 14d71ae5a4SJacob Faibussowitsch PETSC_EXTERN 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) { 50ec1892c8SHong Zhang if (allowzeropivot) { 519566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Zero pivot, row %" PetscInt_FMT "\n", k - 1)); 5270d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE; 5398921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row %" PetscInt_FMT, k - 1); 54ec1892c8SHong Zhang } else { 55943c8ff5SShri Abhyankar a[l + k3] = shift; 56943c8ff5SShri Abhyankar } 57943c8ff5SShri Abhyankar } 5884643e36SBarry Smith 5984643e36SBarry Smith /* interchange if necessary */ 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 stmp = -1. / a[k4]; 6884643e36SBarry Smith i__2 = 2 - k; 6984643e36SBarry Smith aa = &a[1 + k4]; 7026fbe8dcSKarl Rupp for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp; 7184643e36SBarry Smith 7284643e36SBarry Smith /* row elimination with column indexing */ 7384643e36SBarry Smith ax = &a[k4 + 1]; 7484643e36SBarry Smith for (j = kp1; j <= 2; ++j) { 7584643e36SBarry Smith j3 = 2 * j; 7684643e36SBarry Smith stmp = a[l + j3]; 7784643e36SBarry Smith if (l != k) { 7884643e36SBarry Smith a[l + j3] = a[k + j3]; 7984643e36SBarry Smith a[k + j3] = stmp; 8084643e36SBarry Smith } 8184643e36SBarry Smith 8284643e36SBarry Smith i__3 = 2 - k; 8384643e36SBarry Smith ay = &a[1 + k + j3]; 8426fbe8dcSKarl Rupp for (ll = 0; ll < i__3; ll++) ay[ll] += stmp * ax[ll]; 8584643e36SBarry Smith } 8626fbe8dcSKarl Rupp 87b48ee343SBarry Smith ipvt[1] = 2; 882e92ee13SHong Zhang if (a[6] == 0.0) { 89c0aa6a63SJacob Faibussowitsch if (PetscLikely(allowzeropivot)) { 909566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Zero pivot, row 1\n")); 9170d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE; 92c0aa6a63SJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row 1"); 932e92ee13SHong Zhang } 9484643e36SBarry Smith 95c9b7c560SHong Zhang /* Now form the inverse */ 9684643e36SBarry Smith /* compute inverse(u) */ 9784643e36SBarry Smith for (k = 1; k <= 2; ++k) { 9884643e36SBarry Smith k3 = 2 * k; 9984643e36SBarry Smith k4 = k3 + k; 10084643e36SBarry Smith a[k4] = 1.0 / a[k4]; 10184643e36SBarry Smith stmp = -a[k4]; 10284643e36SBarry Smith i__2 = k - 1; 10384643e36SBarry Smith aa = &a[k3 + 1]; 10484643e36SBarry Smith for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp; 10584643e36SBarry Smith kp1 = k + 1; 10684643e36SBarry Smith if (2 < kp1) continue; 10784643e36SBarry Smith ax = aa; 10884643e36SBarry Smith for (j = kp1; j <= 2; ++j) { 10984643e36SBarry Smith j3 = 2 * j; 11084643e36SBarry Smith stmp = a[k + j3]; 11184643e36SBarry Smith a[k + j3] = 0.0; 11284643e36SBarry Smith ay = &a[j3 + 1]; 11326fbe8dcSKarl Rupp for (ll = 0; ll < k; ll++) ay[ll] += stmp * ax[ll]; 11484643e36SBarry Smith } 11584643e36SBarry Smith } 11684643e36SBarry Smith 11784643e36SBarry Smith /* form inverse(u)*inverse(l) */ 11884643e36SBarry Smith k = 1; 11984643e36SBarry Smith k3 = 2 * k; 12084643e36SBarry Smith kp1 = k + 1; 12184643e36SBarry Smith aa = a + k3; 12284643e36SBarry Smith for (i = kp1; i <= 2; ++i) { 123b48ee343SBarry Smith work[i - 1] = aa[i]; 12484643e36SBarry Smith aa[i] = 0.0; 12584643e36SBarry Smith } 12684643e36SBarry Smith for (j = kp1; j <= 2; ++j) { 127b48ee343SBarry Smith stmp = work[j - 1]; 12884643e36SBarry Smith ax = &a[2 * j + 1]; 12984643e36SBarry Smith ay = &a[k3 + 1]; 13084643e36SBarry Smith ay[0] += stmp * ax[0]; 13184643e36SBarry Smith ay[1] += stmp * ax[1]; 13284643e36SBarry Smith } 133b48ee343SBarry Smith l = ipvt[k - 1]; 13484643e36SBarry Smith if (l != k) { 13584643e36SBarry Smith ax = &a[k3 + 1]; 13684643e36SBarry Smith ay = &a[2 * l + 1]; 1379371c9d4SSatish Balay stmp = ax[0]; 1389371c9d4SSatish Balay ax[0] = ay[0]; 1399371c9d4SSatish Balay ay[0] = stmp; 1409371c9d4SSatish Balay stmp = ax[1]; 1419371c9d4SSatish Balay ax[1] = ay[1]; 1429371c9d4SSatish Balay ay[1] = stmp; 14384643e36SBarry Smith } 144*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14584643e36SBarry Smith } 14684643e36SBarry Smith 147ec1892c8SHong Zhang /* gaussian elimination with partial pivoting */ 148d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_9(MatScalar *a, PetscReal shift, PetscBool allowzeropivot, PetscBool *zeropivotdetected) 149d71ae5a4SJacob Faibussowitsch { 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]); 1739371c9d4SSatish Balay if (tmp > max) { 1749371c9d4SSatish Balay max = tmp; 1759371c9d4SSatish Balay l = ll + 1; 1769371c9d4SSatish Balay } 1773dfa136dSBarry Smith } 1783dfa136dSBarry Smith l += k - 1; 1793dfa136dSBarry Smith ipvt[k - 1] = l; 1803dfa136dSBarry Smith 181ec1892c8SHong Zhang if (a[l + k3] == 0.0) { 182ec1892c8SHong Zhang if (shift == 0.0) { 183ec1892c8SHong Zhang if (allowzeropivot) { 1849566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Zero pivot, row %" PetscInt_FMT "\n", k - 1)); 18570d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE; 18698921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row %" PetscInt_FMT, k - 1); 187ec1892c8SHong Zhang } else { 188ec1892c8SHong Zhang a[l + k3] = shift; 189ec1892c8SHong Zhang } 190ec1892c8SHong Zhang } 1913dfa136dSBarry Smith 1923dfa136dSBarry Smith /* interchange if necessary */ 1933dfa136dSBarry Smith if (l != k) { 1943dfa136dSBarry Smith stmp = a[l + k3]; 1953dfa136dSBarry Smith a[l + k3] = a[k4]; 1963dfa136dSBarry Smith a[k4] = stmp; 1973dfa136dSBarry Smith } 1983dfa136dSBarry Smith 1993dfa136dSBarry Smith /* compute multipliers */ 2003dfa136dSBarry Smith stmp = -1. / a[k4]; 2013dfa136dSBarry Smith i__2 = 9 - k; 2023dfa136dSBarry Smith aa = &a[1 + k4]; 20326fbe8dcSKarl Rupp for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp; 2043dfa136dSBarry Smith 2053dfa136dSBarry Smith /* row elimination with column indexing */ 2063dfa136dSBarry Smith ax = &a[k4 + 1]; 2073dfa136dSBarry Smith for (j = kp1; j <= 9; ++j) { 2083dfa136dSBarry Smith j3 = 9 * j; 2093dfa136dSBarry Smith stmp = a[l + j3]; 2103dfa136dSBarry Smith if (l != k) { 2113dfa136dSBarry Smith a[l + j3] = a[k + j3]; 2123dfa136dSBarry Smith a[k + j3] = stmp; 2133dfa136dSBarry Smith } 2143dfa136dSBarry Smith 2153dfa136dSBarry Smith i__3 = 9 - k; 2163dfa136dSBarry Smith ay = &a[1 + k + j3]; 21726fbe8dcSKarl Rupp for (ll = 0; ll < i__3; ll++) ay[ll] += stmp * ax[ll]; 2183dfa136dSBarry Smith } 2193dfa136dSBarry Smith } 2203dfa136dSBarry Smith ipvt[8] = 9; 2212e92ee13SHong Zhang if (a[90] == 0.0) { 222c0aa6a63SJacob Faibussowitsch if (PetscLikely(allowzeropivot)) { 2239566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Zero pivot, row 8\n")); 22470d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE; 225c0aa6a63SJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row 8"); 2262e92ee13SHong Zhang } 2273dfa136dSBarry Smith 228ec1892c8SHong Zhang /* Now form the inverse */ 2293dfa136dSBarry Smith /* compute inverse(u) */ 2303dfa136dSBarry Smith for (k = 1; k <= 9; ++k) { 2313dfa136dSBarry Smith k3 = 9 * k; 2323dfa136dSBarry Smith k4 = k3 + k; 2333dfa136dSBarry Smith a[k4] = 1.0 / a[k4]; 2343dfa136dSBarry Smith stmp = -a[k4]; 2353dfa136dSBarry Smith i__2 = k - 1; 2363dfa136dSBarry Smith aa = &a[k3 + 1]; 2373dfa136dSBarry Smith for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp; 2383dfa136dSBarry Smith kp1 = k + 1; 2393dfa136dSBarry Smith if (9 < kp1) continue; 2403dfa136dSBarry Smith ax = aa; 2413dfa136dSBarry Smith for (j = kp1; j <= 9; ++j) { 2423dfa136dSBarry Smith j3 = 9 * j; 2433dfa136dSBarry Smith stmp = a[k + j3]; 2443dfa136dSBarry Smith a[k + j3] = 0.0; 2453dfa136dSBarry Smith ay = &a[j3 + 1]; 24626fbe8dcSKarl Rupp for (ll = 0; ll < k; ll++) ay[ll] += stmp * ax[ll]; 2473dfa136dSBarry Smith } 2483dfa136dSBarry Smith } 2493dfa136dSBarry Smith 2503dfa136dSBarry Smith /* form inverse(u)*inverse(l) */ 2513dfa136dSBarry Smith for (kb = 1; kb <= 8; ++kb) { 2523dfa136dSBarry Smith k = 9 - kb; 2533dfa136dSBarry Smith k3 = 9 * k; 2543dfa136dSBarry Smith kp1 = k + 1; 2553dfa136dSBarry Smith aa = a + k3; 2563dfa136dSBarry Smith for (i = kp1; i <= 9; ++i) { 2573dfa136dSBarry Smith work[i - 1] = aa[i]; 2583dfa136dSBarry Smith aa[i] = 0.0; 2593dfa136dSBarry Smith } 2603dfa136dSBarry Smith for (j = kp1; j <= 9; ++j) { 2613dfa136dSBarry Smith stmp = work[j - 1]; 2623dfa136dSBarry Smith ax = &a[9 * j + 1]; 2633dfa136dSBarry Smith ay = &a[k3 + 1]; 2643dfa136dSBarry Smith ay[0] += stmp * ax[0]; 2653dfa136dSBarry Smith ay[1] += stmp * ax[1]; 2663dfa136dSBarry Smith ay[2] += stmp * ax[2]; 2673dfa136dSBarry Smith ay[3] += stmp * ax[3]; 2683dfa136dSBarry Smith ay[4] += stmp * ax[4]; 2693dfa136dSBarry Smith ay[5] += stmp * ax[5]; 2703dfa136dSBarry Smith ay[6] += stmp * ax[6]; 2713dfa136dSBarry Smith ay[7] += stmp * ax[7]; 2723dfa136dSBarry Smith ay[8] += stmp * ax[8]; 2733dfa136dSBarry Smith } 2743dfa136dSBarry Smith l = ipvt[k - 1]; 2753dfa136dSBarry Smith if (l != k) { 2763dfa136dSBarry Smith ax = &a[k3 + 1]; 2773dfa136dSBarry Smith ay = &a[9 * l + 1]; 2789371c9d4SSatish Balay stmp = ax[0]; 2799371c9d4SSatish Balay ax[0] = ay[0]; 2809371c9d4SSatish Balay ay[0] = stmp; 2819371c9d4SSatish Balay stmp = ax[1]; 2829371c9d4SSatish Balay ax[1] = ay[1]; 2839371c9d4SSatish Balay ay[1] = stmp; 2849371c9d4SSatish Balay stmp = ax[2]; 2859371c9d4SSatish Balay ax[2] = ay[2]; 2869371c9d4SSatish Balay ay[2] = stmp; 2879371c9d4SSatish Balay stmp = ax[3]; 2889371c9d4SSatish Balay ax[3] = ay[3]; 2899371c9d4SSatish Balay ay[3] = stmp; 2909371c9d4SSatish Balay stmp = ax[4]; 2919371c9d4SSatish Balay ax[4] = ay[4]; 2929371c9d4SSatish Balay ay[4] = stmp; 2939371c9d4SSatish Balay stmp = ax[5]; 2949371c9d4SSatish Balay ax[5] = ay[5]; 2959371c9d4SSatish Balay ay[5] = stmp; 2969371c9d4SSatish Balay stmp = ax[6]; 2979371c9d4SSatish Balay ax[6] = ay[6]; 2989371c9d4SSatish Balay ay[6] = stmp; 2999371c9d4SSatish Balay stmp = ax[7]; 3009371c9d4SSatish Balay ax[7] = ay[7]; 3019371c9d4SSatish Balay ay[7] = stmp; 3029371c9d4SSatish Balay stmp = ax[8]; 3039371c9d4SSatish Balay ax[8] = ay[8]; 3049371c9d4SSatish Balay ay[8] = stmp; 3053dfa136dSBarry Smith } 3063dfa136dSBarry Smith } 307*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3083dfa136dSBarry Smith } 30984643e36SBarry Smith 31029a97285SShri Abhyankar /* 311ec1892c8SHong Zhang Inverts 15 by 15 matrix using gaussian elimination with partial pivoting. 31229a97285SShri Abhyankar 31329a97285SShri Abhyankar Used by the sparse factorization routines in 31429a97285SShri Abhyankar src/mat/impls/baij/seq 31529a97285SShri Abhyankar 31629a97285SShri Abhyankar This is a combination of the Linpack routines 31729a97285SShri Abhyankar dgefa() and dgedi() specialized for a size of 15. 31829a97285SShri Abhyankar 31929a97285SShri Abhyankar */ 32029a97285SShri Abhyankar 321d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_15(MatScalar *a, PetscInt *ipvt, MatScalar *work, PetscReal shift, PetscBool allowzeropivot, PetscBool *zeropivotdetected) 322d71ae5a4SJacob Faibussowitsch { 323766f9fbaSBarry Smith PetscInt i__2, i__3, kp1, j, k, l, ll, i, kb, k3; 32429a97285SShri Abhyankar PetscInt k4, j3; 325766f9fbaSBarry Smith MatScalar *aa, *ax, *ay, stmp; 32629a97285SShri Abhyankar MatReal tmp, max; 32729a97285SShri Abhyankar 32829a97285SShri Abhyankar PetscFunctionBegin; 329c80103daSHong Zhang if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE; 330c80103daSHong Zhang 33129a97285SShri Abhyankar /* Parameter adjustments */ 33229a97285SShri Abhyankar a -= 16; 33329a97285SShri Abhyankar 33429a97285SShri Abhyankar for (k = 1; k <= 14; ++k) { 33529a97285SShri Abhyankar kp1 = k + 1; 33629a97285SShri Abhyankar k3 = 15 * k; 33729a97285SShri Abhyankar k4 = k3 + k; 33829a97285SShri Abhyankar 339ec1892c8SHong Zhang /* find l = pivot index */ 34029a97285SShri Abhyankar i__2 = 16 - k; 34129a97285SShri Abhyankar aa = &a[k4]; 34229a97285SShri Abhyankar max = PetscAbsScalar(aa[0]); 34329a97285SShri Abhyankar l = 1; 34429a97285SShri Abhyankar for (ll = 1; ll < i__2; ll++) { 34529a97285SShri Abhyankar tmp = PetscAbsScalar(aa[ll]); 3469371c9d4SSatish Balay if (tmp > max) { 3479371c9d4SSatish Balay max = tmp; 3489371c9d4SSatish Balay l = ll + 1; 3499371c9d4SSatish Balay } 35029a97285SShri Abhyankar } 35129a97285SShri Abhyankar l += k - 1; 35229a97285SShri Abhyankar ipvt[k - 1] = l; 35329a97285SShri Abhyankar 354ec1892c8SHong Zhang if (a[l + k3] == 0.0) { 355ec1892c8SHong Zhang if (shift == 0.0) { 356ec1892c8SHong Zhang if (allowzeropivot) { 3579566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Zero pivot, row %" PetscInt_FMT "\n", k - 1)); 35870d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE; 35998921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row %" PetscInt_FMT, k - 1); 360ec1892c8SHong Zhang } else { 361ec1892c8SHong Zhang a[l + k3] = shift; 362ec1892c8SHong Zhang } 363ec1892c8SHong Zhang } 36429a97285SShri Abhyankar 36529a97285SShri Abhyankar /* interchange if necessary */ 36629a97285SShri Abhyankar if (l != k) { 36729a97285SShri Abhyankar stmp = a[l + k3]; 36829a97285SShri Abhyankar a[l + k3] = a[k4]; 36929a97285SShri Abhyankar a[k4] = stmp; 37029a97285SShri Abhyankar } 37129a97285SShri Abhyankar 37229a97285SShri Abhyankar /* compute multipliers */ 37329a97285SShri Abhyankar stmp = -1. / a[k4]; 37429a97285SShri Abhyankar i__2 = 15 - k; 37529a97285SShri Abhyankar aa = &a[1 + k4]; 37626fbe8dcSKarl Rupp for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp; 37729a97285SShri Abhyankar 37829a97285SShri Abhyankar /* row elimination with column indexing */ 37929a97285SShri Abhyankar ax = &a[k4 + 1]; 38029a97285SShri Abhyankar for (j = kp1; j <= 15; ++j) { 38129a97285SShri Abhyankar j3 = 15 * j; 38229a97285SShri Abhyankar stmp = a[l + j3]; 38329a97285SShri Abhyankar if (l != k) { 38429a97285SShri Abhyankar a[l + j3] = a[k + j3]; 38529a97285SShri Abhyankar a[k + j3] = stmp; 38629a97285SShri Abhyankar } 38729a97285SShri Abhyankar 38829a97285SShri Abhyankar i__3 = 15 - k; 38929a97285SShri Abhyankar ay = &a[1 + k + j3]; 39026fbe8dcSKarl Rupp for (ll = 0; ll < i__3; ll++) ay[ll] += stmp * ax[ll]; 39129a97285SShri Abhyankar } 39229a97285SShri Abhyankar } 39329a97285SShri Abhyankar ipvt[14] = 15; 394922032d7SHong Zhang if (a[240] == 0.0) { 395c0aa6a63SJacob Faibussowitsch if (PetscLikely(allowzeropivot)) { 3969566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Zero pivot, row 14\n")); 39770d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE; 398c0aa6a63SJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row 14"); 399922032d7SHong Zhang } 40029a97285SShri Abhyankar 401ec1892c8SHong Zhang /* Now form the inverse */ 40229a97285SShri Abhyankar /* compute inverse(u) */ 40329a97285SShri Abhyankar for (k = 1; k <= 15; ++k) { 40429a97285SShri Abhyankar k3 = 15 * k; 40529a97285SShri Abhyankar k4 = k3 + k; 40629a97285SShri Abhyankar a[k4] = 1.0 / a[k4]; 40729a97285SShri Abhyankar stmp = -a[k4]; 40829a97285SShri Abhyankar i__2 = k - 1; 40929a97285SShri Abhyankar aa = &a[k3 + 1]; 41029a97285SShri Abhyankar for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp; 41129a97285SShri Abhyankar kp1 = k + 1; 41229a97285SShri Abhyankar if (15 < kp1) continue; 41329a97285SShri Abhyankar ax = aa; 41429a97285SShri Abhyankar for (j = kp1; j <= 15; ++j) { 41529a97285SShri Abhyankar j3 = 15 * j; 41629a97285SShri Abhyankar stmp = a[k + j3]; 41729a97285SShri Abhyankar a[k + j3] = 0.0; 41829a97285SShri Abhyankar ay = &a[j3 + 1]; 41926fbe8dcSKarl Rupp for (ll = 0; ll < k; ll++) ay[ll] += stmp * ax[ll]; 42029a97285SShri Abhyankar } 42129a97285SShri Abhyankar } 42229a97285SShri Abhyankar 42329a97285SShri Abhyankar /* form inverse(u)*inverse(l) */ 42429a97285SShri Abhyankar for (kb = 1; kb <= 14; ++kb) { 42529a97285SShri Abhyankar k = 15 - kb; 42629a97285SShri Abhyankar k3 = 15 * k; 42729a97285SShri Abhyankar kp1 = k + 1; 42829a97285SShri Abhyankar aa = a + k3; 42929a97285SShri Abhyankar for (i = kp1; i <= 15; ++i) { 43029a97285SShri Abhyankar work[i - 1] = aa[i]; 43129a97285SShri Abhyankar aa[i] = 0.0; 43229a97285SShri Abhyankar } 43329a97285SShri Abhyankar for (j = kp1; j <= 15; ++j) { 43429a97285SShri Abhyankar stmp = work[j - 1]; 43529a97285SShri Abhyankar ax = &a[15 * j + 1]; 43629a97285SShri Abhyankar ay = &a[k3 + 1]; 43729a97285SShri Abhyankar ay[0] += stmp * ax[0]; 43829a97285SShri Abhyankar ay[1] += stmp * ax[1]; 43929a97285SShri Abhyankar ay[2] += stmp * ax[2]; 44029a97285SShri Abhyankar ay[3] += stmp * ax[3]; 44129a97285SShri Abhyankar ay[4] += stmp * ax[4]; 44229a97285SShri Abhyankar ay[5] += stmp * ax[5]; 44329a97285SShri Abhyankar ay[6] += stmp * ax[6]; 44429a97285SShri Abhyankar ay[7] += stmp * ax[7]; 44529a97285SShri Abhyankar ay[8] += stmp * ax[8]; 44629a97285SShri Abhyankar ay[9] += stmp * ax[9]; 44729a97285SShri Abhyankar ay[10] += stmp * ax[10]; 44829a97285SShri Abhyankar ay[11] += stmp * ax[11]; 44929a97285SShri Abhyankar ay[12] += stmp * ax[12]; 45029a97285SShri Abhyankar ay[13] += stmp * ax[13]; 45129a97285SShri Abhyankar ay[14] += stmp * ax[14]; 45229a97285SShri Abhyankar } 45329a97285SShri Abhyankar l = ipvt[k - 1]; 45429a97285SShri Abhyankar if (l != k) { 45529a97285SShri Abhyankar ax = &a[k3 + 1]; 45629a97285SShri Abhyankar ay = &a[15 * l + 1]; 4579371c9d4SSatish Balay stmp = ax[0]; 4589371c9d4SSatish Balay ax[0] = ay[0]; 4599371c9d4SSatish Balay ay[0] = stmp; 4609371c9d4SSatish Balay stmp = ax[1]; 4619371c9d4SSatish Balay ax[1] = ay[1]; 4629371c9d4SSatish Balay ay[1] = stmp; 4639371c9d4SSatish Balay stmp = ax[2]; 4649371c9d4SSatish Balay ax[2] = ay[2]; 4659371c9d4SSatish Balay ay[2] = stmp; 4669371c9d4SSatish Balay stmp = ax[3]; 4679371c9d4SSatish Balay ax[3] = ay[3]; 4689371c9d4SSatish Balay ay[3] = stmp; 4699371c9d4SSatish Balay stmp = ax[4]; 4709371c9d4SSatish Balay ax[4] = ay[4]; 4719371c9d4SSatish Balay ay[4] = stmp; 4729371c9d4SSatish Balay stmp = ax[5]; 4739371c9d4SSatish Balay ax[5] = ay[5]; 4749371c9d4SSatish Balay ay[5] = stmp; 4759371c9d4SSatish Balay stmp = ax[6]; 4769371c9d4SSatish Balay ax[6] = ay[6]; 4779371c9d4SSatish Balay ay[6] = stmp; 4789371c9d4SSatish Balay stmp = ax[7]; 4799371c9d4SSatish Balay ax[7] = ay[7]; 4809371c9d4SSatish Balay ay[7] = stmp; 4819371c9d4SSatish Balay stmp = ax[8]; 4829371c9d4SSatish Balay ax[8] = ay[8]; 4839371c9d4SSatish Balay ay[8] = stmp; 4849371c9d4SSatish Balay stmp = ax[9]; 4859371c9d4SSatish Balay ax[9] = ay[9]; 4869371c9d4SSatish Balay ay[9] = stmp; 4879371c9d4SSatish Balay stmp = ax[10]; 4889371c9d4SSatish Balay ax[10] = ay[10]; 4899371c9d4SSatish Balay ay[10] = stmp; 4909371c9d4SSatish Balay stmp = ax[11]; 4919371c9d4SSatish Balay ax[11] = ay[11]; 4929371c9d4SSatish Balay ay[11] = stmp; 4939371c9d4SSatish Balay stmp = ax[12]; 4949371c9d4SSatish Balay ax[12] = ay[12]; 4959371c9d4SSatish Balay ay[12] = stmp; 4969371c9d4SSatish Balay stmp = ax[13]; 4979371c9d4SSatish Balay ax[13] = ay[13]; 4989371c9d4SSatish Balay ay[13] = stmp; 4999371c9d4SSatish Balay stmp = ax[14]; 5009371c9d4SSatish Balay ax[14] = ay[14]; 5019371c9d4SSatish Balay ay[14] = stmp; 50229a97285SShri Abhyankar } 50329a97285SShri Abhyankar } 504*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 50529a97285SShri Abhyankar } 506