1be1d678aSKris Buschelman 2ede5e988SBarry Smith /* 3ec1892c8SHong Zhang Inverts 3 by 3 matrix using gaussian elimination with partial pivoting. 471c5468dSBarry Smith 571c5468dSBarry Smith Used by the sparse factorization routines in 6dd882469SBarry Smith src/mat/impls/baij/seq 771c5468dSBarry Smith 871c5468dSBarry Smith This is a combination of the Linpack routines 971c5468dSBarry Smith dgefa() and dgedi() specialized for a size of 3. 1071c5468dSBarry Smith 11ede5e988SBarry Smith */ 12c6db04a5SJed Brown #include <petscsys.h> 13ede5e988SBarry Smith 14*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_3(MatScalar *a, PetscReal shift, PetscBool allowzeropivot, PetscBool *zeropivotdetected) 15*d71ae5a4SJacob Faibussowitsch { 16690b6cddSBarry Smith PetscInt i__2, i__3, kp1, j, k, l, ll, i, ipvt[3], kb, k3; 17690b6cddSBarry Smith PetscInt k4, j3; 18b48ee343SBarry Smith MatScalar *aa, *ax, *ay, work[9], stmp; 19329f5518SBarry Smith MatReal tmp, max; 20ede5e988SBarry Smith 213a40ed3dSBarry Smith PetscFunctionBegin; 22c80103daSHong Zhang if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE; 230daa89e9SHong Zhang shift = .333 * shift * (1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[4]) + PetscAbsScalar(a[8])); 24ec1892c8SHong Zhang 25ede5e988SBarry Smith /* Parameter adjustments */ 26ede5e988SBarry Smith a -= 4; 27ede5e988SBarry Smith 28ede5e988SBarry Smith for (k = 1; k <= 2; ++k) { 29ede5e988SBarry Smith kp1 = k + 1; 3073d4a2d6SBarry Smith k3 = 3 * k; 3173d4a2d6SBarry Smith k4 = k3 + k; 32ede5e988SBarry Smith 33ec1892c8SHong Zhang /* find l = pivot index */ 34ede5e988SBarry Smith i__2 = 4 - k; 3573d4a2d6SBarry Smith aa = &a[k4]; 36ede5e988SBarry Smith max = PetscAbsScalar(aa[0]); 37ede5e988SBarry Smith l = 1; 38ede5e988SBarry Smith for (ll = 1; ll < i__2; ll++) { 39ede5e988SBarry Smith tmp = PetscAbsScalar(aa[ll]); 409371c9d4SSatish Balay if (tmp > max) { 419371c9d4SSatish Balay max = tmp; 429371c9d4SSatish Balay l = ll + 1; 439371c9d4SSatish Balay } 44ede5e988SBarry Smith } 45ede5e988SBarry Smith l += k - 1; 46da10e913SBarry Smith ipvt[k - 1] = l; 47ede5e988SBarry Smith 48ab055debSShri 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 { 55ab055debSShri Abhyankar /* Shift is applied to single diagonal entry */ 56ab055debSShri Abhyankar a[l + k3] = shift; 57ab055debSShri Abhyankar } 58ab055debSShri Abhyankar } 59ede5e988SBarry Smith 60ec1892c8SHong Zhang /* interchange if necessary */ 61ede5e988SBarry Smith if (l != k) { 6239d66777SBarry Smith stmp = a[l + k3]; 6373d4a2d6SBarry Smith a[l + k3] = a[k4]; 6439d66777SBarry Smith a[k4] = stmp; 65ede5e988SBarry Smith } 66ede5e988SBarry Smith 67ede5e988SBarry Smith /* compute multipliers */ 6839d66777SBarry Smith stmp = -1. / a[k4]; 69ede5e988SBarry Smith i__2 = 3 - k; 7073d4a2d6SBarry Smith aa = &a[1 + k4]; 7126fbe8dcSKarl Rupp for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp; 72ede5e988SBarry Smith 73ede5e988SBarry Smith /* row elimination with column indexing */ 7473d4a2d6SBarry Smith ax = &a[k4 + 1]; 75ede5e988SBarry Smith for (j = kp1; j <= 3; ++j) { 7673d4a2d6SBarry Smith j3 = 3 * j; 7739d66777SBarry Smith stmp = a[l + j3]; 78ede5e988SBarry Smith if (l != k) { 7973d4a2d6SBarry Smith a[l + j3] = a[k + j3]; 8039d66777SBarry Smith a[k + j3] = stmp; 81ede5e988SBarry Smith } 82ede5e988SBarry Smith 83ede5e988SBarry Smith i__3 = 3 - k; 8439d66777SBarry Smith ay = &a[1 + k + j3]; 8526fbe8dcSKarl Rupp for (ll = 0; ll < i__3; ll++) ay[ll] += stmp * ax[ll]; 86ede5e988SBarry Smith } 87ede5e988SBarry Smith } 88da10e913SBarry Smith ipvt[2] = 3; 89cd545bacSHong Zhang if (a[12] == 0.0) { 90c0aa6a63SJacob Faibussowitsch if (PetscLikely(allowzeropivot)) { 919566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Zero pivot, row 2\n")); 9270d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE; 93c0aa6a63SJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row 2"); 94cd545bacSHong Zhang } 95ede5e988SBarry Smith 96ec1892c8SHong Zhang /* Now form the inverse */ 97ede5e988SBarry Smith /* compute inverse(u) */ 9825783f72SBarry Smith for (k = 1; k <= 3; ++k) { 9973d4a2d6SBarry Smith k3 = 3 * k; 10073d4a2d6SBarry Smith k4 = k3 + k; 10173d4a2d6SBarry Smith a[k4] = 1.0 / a[k4]; 10239d66777SBarry Smith stmp = -a[k4]; 103ede5e988SBarry Smith i__2 = k - 1; 10473d4a2d6SBarry Smith aa = &a[k3 + 1]; 10539d66777SBarry Smith for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp; 106ede5e988SBarry Smith kp1 = k + 1; 10725783f72SBarry Smith if (3 < kp1) continue; 108ede5e988SBarry Smith ax = aa; 10925783f72SBarry Smith for (j = kp1; j <= 3; ++j) { 11073d4a2d6SBarry Smith j3 = 3 * j; 11139d66777SBarry Smith stmp = a[k + j3]; 11239d66777SBarry Smith a[k + j3] = 0.0; 11373d4a2d6SBarry Smith ay = &a[j3 + 1]; 11426fbe8dcSKarl Rupp for (ll = 0; ll < k; ll++) ay[ll] += stmp * ax[ll]; 115ede5e988SBarry Smith } 116ede5e988SBarry Smith } 117ede5e988SBarry Smith 118ede5e988SBarry Smith /* form inverse(u)*inverse(l) */ 11925783f72SBarry Smith for (kb = 1; kb <= 2; ++kb) { 12025783f72SBarry Smith k = 3 - kb; 12173d4a2d6SBarry Smith k3 = 3 * k; 122ede5e988SBarry Smith kp1 = k + 1; 12373d4a2d6SBarry Smith aa = a + k3; 12425783f72SBarry Smith for (i = kp1; i <= 3; ++i) { 125b48ee343SBarry Smith work[i - 1] = aa[i]; 12673d4a2d6SBarry Smith aa[i] = 0.0; 127ede5e988SBarry Smith } 12825783f72SBarry Smith for (j = kp1; j <= 3; ++j) { 129b48ee343SBarry Smith stmp = work[j - 1]; 13025783f72SBarry Smith ax = &a[3 * j + 1]; 13173d4a2d6SBarry Smith ay = &a[k3 + 1]; 13239d66777SBarry Smith ay[0] += stmp * ax[0]; 13339d66777SBarry Smith ay[1] += stmp * ax[1]; 13439d66777SBarry Smith ay[2] += stmp * ax[2]; 135ede5e988SBarry Smith } 136da10e913SBarry Smith l = ipvt[k - 1]; 137ede5e988SBarry Smith if (l != k) { 13873d4a2d6SBarry Smith ax = &a[k3 + 1]; 13925783f72SBarry Smith ay = &a[3 * 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; 149ede5e988SBarry Smith } 150ede5e988SBarry Smith } 1513a40ed3dSBarry Smith PetscFunctionReturn(0); 152ede5e988SBarry Smith } 153