1be1d678aSKris Buschelman 284643e36SBarry Smith /* 3ec1892c8SHong Zhang Inverts 7 by 7 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 7. 1084643e36SBarry Smith 1184643e36SBarry Smith */ 12c6db04a5SJed Brown #include <petscsys.h> 1384643e36SBarry Smith 14*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_7(MatScalar *a, PetscReal shift, PetscBool allowzeropivot, PetscBool *zeropivotdetected) { 15690b6cddSBarry Smith PetscInt i__2, i__3, kp1, j, k, l, ll, i, ipvt[7], kb, k3; 16690b6cddSBarry Smith PetscInt k4, j3; 171f1046a5SBarry Smith MatScalar *aa, *ax, *ay, work[49], 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[8]) + PetscAbsScalar(a[16]) + PetscAbsScalar(a[24]) + PetscAbsScalar(a[32]) + PetscAbsScalar(a[40]) + PetscAbsScalar(a[48])); 23943c8ff5SShri Abhyankar 2484643e36SBarry Smith /* Parameter adjustments */ 2584643e36SBarry Smith a -= 8; 2684643e36SBarry Smith 2784643e36SBarry Smith for (k = 1; k <= 6; ++k) { 2884643e36SBarry Smith kp1 = k + 1; 2984643e36SBarry Smith k3 = 7 * k; 3084643e36SBarry Smith k4 = k3 + k; 3184643e36SBarry Smith 32ec1892c8SHong Zhang /* find l = pivot index */ 33ed33f8a5SSatish Balay i__2 = 8 - 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; 451f1046a5SBarry 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 /* SHIFT is applied to SINGLE diagonal entry; does this make any sense? */ 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 = 7 - 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 <= 7; ++j) { 7584643e36SBarry Smith j3 = 7 * 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 = 7 - k; 8384643e36SBarry Smith ay = &a[1 + k + j3]; 8426fbe8dcSKarl Rupp for (ll = 0; ll < i__3; ll++) ay[ll] += stmp * ax[ll]; 8584643e36SBarry Smith } 8684643e36SBarry Smith } 871f1046a5SBarry Smith ipvt[6] = 7; 882e92ee13SHong Zhang if (a[56] == 0.0) { 89c0aa6a63SJacob Faibussowitsch if (PetscLikely(allowzeropivot)) { 909566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Zero pivot, row 6\n")); 9170d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE; 92c0aa6a63SJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row 6"); 932e92ee13SHong Zhang } 9484643e36SBarry Smith 95ec1892c8SHong Zhang /* Now form the inverse */ 9684643e36SBarry Smith /* compute inverse(u) */ 9784643e36SBarry Smith for (k = 1; k <= 7; ++k) { 9884643e36SBarry Smith k3 = 7 * 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 (7 < kp1) continue; 10784643e36SBarry Smith ax = aa; 10884643e36SBarry Smith for (j = kp1; j <= 7; ++j) { 10984643e36SBarry Smith j3 = 7 * 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 for (kb = 1; kb <= 6; ++kb) { 11984643e36SBarry Smith k = 7 - kb; 12084643e36SBarry Smith k3 = 7 * k; 12184643e36SBarry Smith kp1 = k + 1; 12284643e36SBarry Smith aa = a + k3; 12384643e36SBarry Smith for (i = kp1; i <= 7; ++i) { 1241f1046a5SBarry Smith work[i - 1] = aa[i]; 12584643e36SBarry Smith aa[i] = 0.0; 12684643e36SBarry Smith } 12784643e36SBarry Smith for (j = kp1; j <= 7; ++j) { 1281f1046a5SBarry Smith stmp = work[j - 1]; 12984643e36SBarry Smith ax = &a[7 * j + 1]; 13084643e36SBarry Smith ay = &a[k3 + 1]; 13184643e36SBarry Smith ay[0] += stmp * ax[0]; 13284643e36SBarry Smith ay[1] += stmp * ax[1]; 13384643e36SBarry Smith ay[2] += stmp * ax[2]; 13484643e36SBarry Smith ay[3] += stmp * ax[3]; 13584643e36SBarry Smith ay[4] += stmp * ax[4]; 13684643e36SBarry Smith ay[5] += stmp * ax[5]; 13784643e36SBarry Smith ay[6] += stmp * ax[6]; 13884643e36SBarry Smith } 1391f1046a5SBarry Smith l = ipvt[k - 1]; 14084643e36SBarry Smith if (l != k) { 14184643e36SBarry Smith ax = &a[k3 + 1]; 14215091d37SBarry Smith ay = &a[7 * l + 1]; 143*9371c9d4SSatish Balay stmp = ax[0]; 144*9371c9d4SSatish Balay ax[0] = ay[0]; 145*9371c9d4SSatish Balay ay[0] = stmp; 146*9371c9d4SSatish Balay stmp = ax[1]; 147*9371c9d4SSatish Balay ax[1] = ay[1]; 148*9371c9d4SSatish Balay ay[1] = stmp; 149*9371c9d4SSatish Balay stmp = ax[2]; 150*9371c9d4SSatish Balay ax[2] = ay[2]; 151*9371c9d4SSatish Balay ay[2] = stmp; 152*9371c9d4SSatish Balay stmp = ax[3]; 153*9371c9d4SSatish Balay ax[3] = ay[3]; 154*9371c9d4SSatish Balay ay[3] = stmp; 155*9371c9d4SSatish Balay stmp = ax[4]; 156*9371c9d4SSatish Balay ax[4] = ay[4]; 157*9371c9d4SSatish Balay ay[4] = stmp; 158*9371c9d4SSatish Balay stmp = ax[5]; 159*9371c9d4SSatish Balay ax[5] = ay[5]; 160*9371c9d4SSatish Balay ay[5] = stmp; 161*9371c9d4SSatish Balay stmp = ax[6]; 162*9371c9d4SSatish Balay ax[6] = ay[6]; 163*9371c9d4SSatish Balay ay[6] = stmp; 16484643e36SBarry Smith } 16584643e36SBarry Smith } 16684643e36SBarry Smith PetscFunctionReturn(0); 16784643e36SBarry Smith } 168