14224c193SBarry Smith /* 2ec1892c8SHong Zhang Inverts 5 by 5 matrix using gaussian elimination with partial pivoting. 371c5468dSBarry Smith 471c5468dSBarry Smith Used by the sparse factorization routines in 5dd882469SBarry Smith src/mat/impls/baij/seq 671c5468dSBarry Smith 771c5468dSBarry Smith This is a combination of the Linpack routines 871c5468dSBarry Smith dgefa() and dgedi() specialized for a size of 5. 971c5468dSBarry Smith 104224c193SBarry Smith */ 11c6db04a5SJed Brown #include <petscsys.h> 12d6acfc2dSPierre Jolivet #include <petsc/private/kernels/blockinvert.h> 134224c193SBarry Smith 14d6acfc2dSPierre Jolivet PetscErrorCode PetscKernel_A_gets_inverse_A_5(MatScalar *a, PetscInt *ipvt, MatScalar *work, PetscReal shift, PetscBool allowzeropivot, PetscBool *zeropivotdetected) 15d71ae5a4SJacob Faibussowitsch { 16766f9fbaSBarry Smith PetscInt i__2, i__3, kp1, j, k, l, ll, i, kb, k3; 17690b6cddSBarry Smith PetscInt k4, j3; 18766f9fbaSBarry Smith MatScalar *aa, *ax, *ay, stmp; 19329f5518SBarry Smith MatReal tmp, max; 204224c193SBarry Smith 213a40ed3dSBarry Smith PetscFunctionBegin; 22c80103daSHong Zhang if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE; 23943c8ff5SShri Abhyankar shift = .25 * shift * (1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[6]) + PetscAbsScalar(a[12]) + PetscAbsScalar(a[18]) + PetscAbsScalar(a[24])); 24ec1892c8SHong Zhang 254224c193SBarry Smith /* Parameter adjustments */ 268a36c062SBarry Smith a -= 6; 274224c193SBarry Smith 288a36c062SBarry Smith for (k = 1; k <= 4; ++k) { 294224c193SBarry Smith kp1 = k + 1; 308a36c062SBarry Smith k3 = 5 * k; 314224c193SBarry Smith k4 = k3 + k; 324224c193SBarry Smith 33ec1892c8SHong Zhang /* find l = pivot index */ 34ed33f8a5SSatish Balay i__2 = 6 - k; 354224c193SBarry Smith aa = &a[k4]; 364224c193SBarry Smith max = PetscAbsScalar(aa[0]); 374224c193SBarry Smith l = 1; 384224c193SBarry Smith for (ll = 1; ll < i__2; ll++) { 394224c193SBarry Smith tmp = PetscAbsScalar(aa[ll]); 409371c9d4SSatish Balay if (tmp > max) { 419371c9d4SSatish Balay max = tmp; 429371c9d4SSatish Balay l = ll + 1; 439371c9d4SSatish Balay } 444224c193SBarry Smith } 454224c193SBarry Smith l += k - 1; 46b48ee343SBarry Smith ipvt[k - 1] = l; 474224c193SBarry 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)); 52ec1892c8SHong Zhang *zeropivotdetected = PETSC_TRUE; 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 } 584224c193SBarry Smith 594224c193SBarry Smith /* interchange if necessary */ 604224c193SBarry Smith if (l != k) { 614224c193SBarry Smith stmp = a[l + k3]; 624224c193SBarry Smith a[l + k3] = a[k4]; 634224c193SBarry Smith a[k4] = stmp; 644224c193SBarry Smith } 654224c193SBarry Smith 664224c193SBarry Smith /* compute multipliers */ 674224c193SBarry Smith stmp = -1. / a[k4]; 688a36c062SBarry Smith i__2 = 5 - k; 694224c193SBarry Smith aa = &a[1 + k4]; 7026fbe8dcSKarl Rupp for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp; 714224c193SBarry Smith 724224c193SBarry Smith /* row elimination with column indexing */ 734224c193SBarry Smith ax = &a[k4 + 1]; 748a36c062SBarry Smith for (j = kp1; j <= 5; ++j) { 758a36c062SBarry Smith j3 = 5 * j; 764224c193SBarry Smith stmp = a[l + j3]; 774224c193SBarry Smith if (l != k) { 784224c193SBarry Smith a[l + j3] = a[k + j3]; 794224c193SBarry Smith a[k + j3] = stmp; 804224c193SBarry Smith } 814224c193SBarry Smith 828a36c062SBarry Smith i__3 = 5 - k; 834224c193SBarry Smith ay = &a[1 + k + j3]; 8426fbe8dcSKarl Rupp for (ll = 0; ll < i__3; ll++) ay[ll] += stmp * ax[ll]; 854224c193SBarry Smith } 864224c193SBarry Smith } 87b48ee343SBarry Smith ipvt[4] = 5; 882e92ee13SHong Zhang if (a[30] == 0.0) { 89*966bd95aSPierre Jolivet PetscCheck(PetscLikely(allowzeropivot), PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row 4"); 909566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Zero pivot, row 4\n")); 912e92ee13SHong Zhang *zeropivotdetected = PETSC_TRUE; 922e92ee13SHong Zhang } 934224c193SBarry Smith 94ec1892c8SHong Zhang /* Now form the inverse */ 954224c193SBarry Smith /* compute inverse(u) */ 968a36c062SBarry Smith for (k = 1; k <= 5; ++k) { 978a36c062SBarry Smith k3 = 5 * k; 984224c193SBarry Smith k4 = k3 + k; 994224c193SBarry Smith a[k4] = 1.0 / a[k4]; 1004224c193SBarry Smith stmp = -a[k4]; 1014224c193SBarry Smith i__2 = k - 1; 1024224c193SBarry Smith aa = &a[k3 + 1]; 1034224c193SBarry Smith for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp; 1044224c193SBarry Smith kp1 = k + 1; 1058a36c062SBarry Smith if (5 < kp1) continue; 1064224c193SBarry Smith ax = aa; 1078a36c062SBarry Smith for (j = kp1; j <= 5; ++j) { 1088a36c062SBarry Smith j3 = 5 * j; 1094224c193SBarry Smith stmp = a[k + j3]; 1104224c193SBarry Smith a[k + j3] = 0.0; 1114224c193SBarry Smith ay = &a[j3 + 1]; 11226fbe8dcSKarl Rupp for (ll = 0; ll < k; ll++) ay[ll] += stmp * ax[ll]; 1134224c193SBarry Smith } 1144224c193SBarry Smith } 1154224c193SBarry Smith 1164224c193SBarry Smith /* form inverse(u)*inverse(l) */ 1178a36c062SBarry Smith for (kb = 1; kb <= 4; ++kb) { 1188a36c062SBarry Smith k = 5 - kb; 1198a36c062SBarry Smith k3 = 5 * k; 1204224c193SBarry Smith kp1 = k + 1; 1214224c193SBarry Smith aa = a + k3; 1228a36c062SBarry Smith for (i = kp1; i <= 5; ++i) { 123b48ee343SBarry Smith work[i - 1] = aa[i]; 1244224c193SBarry Smith aa[i] = 0.0; 1254224c193SBarry Smith } 1268a36c062SBarry Smith for (j = kp1; j <= 5; ++j) { 127b48ee343SBarry Smith stmp = work[j - 1]; 1288a36c062SBarry Smith ax = &a[5 * j + 1]; 1294224c193SBarry Smith ay = &a[k3 + 1]; 1304224c193SBarry Smith ay[0] += stmp * ax[0]; 1314224c193SBarry Smith ay[1] += stmp * ax[1]; 1324224c193SBarry Smith ay[2] += stmp * ax[2]; 1338a36c062SBarry Smith ay[3] += stmp * ax[3]; 1348a36c062SBarry Smith ay[4] += stmp * ax[4]; 1354224c193SBarry Smith } 136b48ee343SBarry Smith l = ipvt[k - 1]; 1374224c193SBarry Smith if (l != k) { 1384224c193SBarry Smith ax = &a[k3 + 1]; 1398a36c062SBarry Smith ay = &a[5 * l + 1]; 1409371c9d4SSatish Balay stmp = ax[0]; 1419371c9d4SSatish Balay ax[0] = ay[0]; 1429371c9d4SSatish Balay ay[0] = stmp; 1439371c9d4SSatish Balay stmp = ax[1]; 1449371c9d4SSatish Balay ax[1] = ay[1]; 1459371c9d4SSatish Balay ay[1] = stmp; 1469371c9d4SSatish Balay stmp = ax[2]; 1479371c9d4SSatish Balay ax[2] = ay[2]; 1489371c9d4SSatish Balay ay[2] = stmp; 1499371c9d4SSatish Balay stmp = ax[3]; 1509371c9d4SSatish Balay ax[3] = ay[3]; 1519371c9d4SSatish Balay ay[3] = stmp; 1529371c9d4SSatish Balay stmp = ax[4]; 1539371c9d4SSatish Balay ax[4] = ay[4]; 1549371c9d4SSatish Balay ay[4] = stmp; 1554224c193SBarry Smith } 1564224c193SBarry Smith } 1573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1584224c193SBarry Smith } 159