1ede5e988SBarry Smith /* 2ec1892c8SHong Zhang Inverts 3 by 3 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 3. 971c5468dSBarry Smith 10ede5e988SBarry Smith */ 11c6db04a5SJed Brown #include <petscsys.h> 12d6acfc2dSPierre Jolivet #include <petsc/private/kernels/blockinvert.h> 13ede5e988SBarry Smith 14d6acfc2dSPierre Jolivet PetscErrorCode PetscKernel_A_gets_inverse_A_3(MatScalar *a, PetscReal shift, PetscBool allowzeropivot, PetscBool *zeropivotdetected) 15d71ae5a4SJacob 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) { 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 { 54ab055debSShri Abhyankar /* Shift is applied to single diagonal entry */ 55ab055debSShri Abhyankar a[l + k3] = shift; 56ab055debSShri Abhyankar } 57ab055debSShri Abhyankar } 58ede5e988SBarry Smith 59ec1892c8SHong Zhang /* interchange if necessary */ 60ede5e988SBarry Smith if (l != k) { 6139d66777SBarry Smith stmp = a[l + k3]; 6273d4a2d6SBarry Smith a[l + k3] = a[k4]; 6339d66777SBarry Smith a[k4] = stmp; 64ede5e988SBarry Smith } 65ede5e988SBarry Smith 66ede5e988SBarry Smith /* compute multipliers */ 6739d66777SBarry Smith stmp = -1. / a[k4]; 68ede5e988SBarry Smith i__2 = 3 - k; 6973d4a2d6SBarry Smith aa = &a[1 + k4]; 7026fbe8dcSKarl Rupp for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp; 71ede5e988SBarry Smith 72ede5e988SBarry Smith /* row elimination with column indexing */ 7373d4a2d6SBarry Smith ax = &a[k4 + 1]; 74ede5e988SBarry Smith for (j = kp1; j <= 3; ++j) { 7573d4a2d6SBarry Smith j3 = 3 * j; 7639d66777SBarry Smith stmp = a[l + j3]; 77ede5e988SBarry Smith if (l != k) { 7873d4a2d6SBarry Smith a[l + j3] = a[k + j3]; 7939d66777SBarry Smith a[k + j3] = stmp; 80ede5e988SBarry Smith } 81ede5e988SBarry Smith 82ede5e988SBarry Smith i__3 = 3 - k; 8339d66777SBarry Smith ay = &a[1 + k + j3]; 8426fbe8dcSKarl Rupp for (ll = 0; ll < i__3; ll++) ay[ll] += stmp * ax[ll]; 85ede5e988SBarry Smith } 86ede5e988SBarry Smith } 87da10e913SBarry Smith ipvt[2] = 3; 88cd545bacSHong Zhang if (a[12] == 0.0) { 89*966bd95aSPierre Jolivet PetscCheck(PetscLikely(allowzeropivot), PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row 2"); 909566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Zero pivot, row 2\n")); 9170d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE; 92cd545bacSHong Zhang } 93ede5e988SBarry Smith 94ec1892c8SHong Zhang /* Now form the inverse */ 95ede5e988SBarry Smith /* compute inverse(u) */ 9625783f72SBarry Smith for (k = 1; k <= 3; ++k) { 9773d4a2d6SBarry Smith k3 = 3 * k; 9873d4a2d6SBarry Smith k4 = k3 + k; 9973d4a2d6SBarry Smith a[k4] = 1.0 / a[k4]; 10039d66777SBarry Smith stmp = -a[k4]; 101ede5e988SBarry Smith i__2 = k - 1; 10273d4a2d6SBarry Smith aa = &a[k3 + 1]; 10339d66777SBarry Smith for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp; 104ede5e988SBarry Smith kp1 = k + 1; 10525783f72SBarry Smith if (3 < kp1) continue; 106ede5e988SBarry Smith ax = aa; 10725783f72SBarry Smith for (j = kp1; j <= 3; ++j) { 10873d4a2d6SBarry Smith j3 = 3 * j; 10939d66777SBarry Smith stmp = a[k + j3]; 11039d66777SBarry Smith a[k + j3] = 0.0; 11173d4a2d6SBarry Smith ay = &a[j3 + 1]; 11226fbe8dcSKarl Rupp for (ll = 0; ll < k; ll++) ay[ll] += stmp * ax[ll]; 113ede5e988SBarry Smith } 114ede5e988SBarry Smith } 115ede5e988SBarry Smith 116ede5e988SBarry Smith /* form inverse(u)*inverse(l) */ 11725783f72SBarry Smith for (kb = 1; kb <= 2; ++kb) { 11825783f72SBarry Smith k = 3 - kb; 11973d4a2d6SBarry Smith k3 = 3 * k; 120ede5e988SBarry Smith kp1 = k + 1; 12173d4a2d6SBarry Smith aa = a + k3; 12225783f72SBarry Smith for (i = kp1; i <= 3; ++i) { 123b48ee343SBarry Smith work[i - 1] = aa[i]; 12473d4a2d6SBarry Smith aa[i] = 0.0; 125ede5e988SBarry Smith } 12625783f72SBarry Smith for (j = kp1; j <= 3; ++j) { 127b48ee343SBarry Smith stmp = work[j - 1]; 12825783f72SBarry Smith ax = &a[3 * j + 1]; 12973d4a2d6SBarry Smith ay = &a[k3 + 1]; 13039d66777SBarry Smith ay[0] += stmp * ax[0]; 13139d66777SBarry Smith ay[1] += stmp * ax[1]; 13239d66777SBarry Smith ay[2] += stmp * ax[2]; 133ede5e988SBarry Smith } 134da10e913SBarry Smith l = ipvt[k - 1]; 135ede5e988SBarry Smith if (l != k) { 13673d4a2d6SBarry Smith ax = &a[k3 + 1]; 13725783f72SBarry Smith ay = &a[3 * l + 1]; 1389371c9d4SSatish Balay stmp = ax[0]; 1399371c9d4SSatish Balay ax[0] = ay[0]; 1409371c9d4SSatish Balay ay[0] = stmp; 1419371c9d4SSatish Balay stmp = ax[1]; 1429371c9d4SSatish Balay ax[1] = ay[1]; 1439371c9d4SSatish Balay ay[1] = stmp; 1449371c9d4SSatish Balay stmp = ax[2]; 1459371c9d4SSatish Balay ax[2] = ay[2]; 1469371c9d4SSatish Balay ay[2] = stmp; 147ede5e988SBarry Smith } 148ede5e988SBarry Smith } 1493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 150ede5e988SBarry Smith } 151