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 146baedc03SHong Zhang PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_3(MatScalar *a,PetscReal shift,PetscBool allowzeropivot,PetscBool *zeropivotdetected) 15ede5e988SBarry Smith { 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]); 40ede5e988SBarry Smith if (tmp > max) { max = tmp; l = ll+1;} 41ede5e988SBarry Smith } 42ede5e988SBarry Smith l += k - 1; 43da10e913SBarry Smith ipvt[k-1] = l; 44ede5e988SBarry Smith 45ab055debSShri Abhyankar if (a[l + k3] == 0.0) { 46ec1892c8SHong Zhang if (shift == 0.0) { 47ec1892c8SHong Zhang if (allowzeropivot) { 48*9566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL,"Zero pivot, row %" PetscInt_FMT "\n",k-1)); 4970d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE; 5098921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %" PetscInt_FMT,k-1); 51ec1892c8SHong Zhang } else { 52ab055debSShri Abhyankar /* Shift is applied to single diagonal entry */ 53ab055debSShri Abhyankar a[l + k3] = shift; 54ab055debSShri Abhyankar } 55ab055debSShri Abhyankar } 56ede5e988SBarry Smith 57ec1892c8SHong Zhang /* interchange if necessary */ 58ede5e988SBarry Smith if (l != k) { 5939d66777SBarry Smith stmp = a[l + k3]; 6073d4a2d6SBarry Smith a[l + k3] = a[k4]; 6139d66777SBarry Smith a[k4] = stmp; 62ede5e988SBarry Smith } 63ede5e988SBarry Smith 64ede5e988SBarry Smith /* compute multipliers */ 6539d66777SBarry Smith stmp = -1. / a[k4]; 66ede5e988SBarry Smith i__2 = 3 - k; 6773d4a2d6SBarry Smith aa = &a[1 + k4]; 6826fbe8dcSKarl Rupp for (ll=0; ll<i__2; ll++) aa[ll] *= stmp; 69ede5e988SBarry Smith 70ede5e988SBarry Smith /* row elimination with column indexing */ 7173d4a2d6SBarry Smith ax = &a[k4+1]; 72ede5e988SBarry Smith for (j = kp1; j <= 3; ++j) { 7373d4a2d6SBarry Smith j3 = 3*j; 7439d66777SBarry Smith stmp = a[l + j3]; 75ede5e988SBarry Smith if (l != k) { 7673d4a2d6SBarry Smith a[l + j3] = a[k + j3]; 7739d66777SBarry Smith a[k + j3] = stmp; 78ede5e988SBarry Smith } 79ede5e988SBarry Smith 80ede5e988SBarry Smith i__3 = 3 - k; 8139d66777SBarry Smith ay = &a[1+k+j3]; 8226fbe8dcSKarl Rupp for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll]; 83ede5e988SBarry Smith } 84ede5e988SBarry Smith } 85da10e913SBarry Smith ipvt[2] = 3; 86cd545bacSHong Zhang if (a[12] == 0.0) { 87c0aa6a63SJacob Faibussowitsch if (PetscLikely(allowzeropivot)) { 88*9566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL,"Zero pivot, row 2\n")); 8970d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE; 90c0aa6a63SJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row 2"); 91cd545bacSHong Zhang } 92ede5e988SBarry Smith 93ec1892c8SHong Zhang /* Now form the inverse */ 94ede5e988SBarry Smith /* compute inverse(u) */ 9525783f72SBarry Smith for (k = 1; k <= 3; ++k) { 9673d4a2d6SBarry Smith k3 = 3*k; 9773d4a2d6SBarry Smith k4 = k3 + k; 9873d4a2d6SBarry Smith a[k4] = 1.0 / a[k4]; 9939d66777SBarry Smith stmp = -a[k4]; 100ede5e988SBarry Smith i__2 = k - 1; 10173d4a2d6SBarry Smith aa = &a[k3 + 1]; 10239d66777SBarry Smith for (ll=0; ll<i__2; ll++) aa[ll] *= stmp; 103ede5e988SBarry Smith kp1 = k + 1; 10425783f72SBarry Smith if (3 < kp1) continue; 105ede5e988SBarry Smith ax = aa; 10625783f72SBarry Smith for (j = kp1; j <= 3; ++j) { 10773d4a2d6SBarry Smith j3 = 3*j; 10839d66777SBarry Smith stmp = a[k + j3]; 10939d66777SBarry Smith a[k + j3] = 0.0; 11073d4a2d6SBarry Smith ay = &a[j3 + 1]; 11126fbe8dcSKarl Rupp for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll]; 112ede5e988SBarry Smith } 113ede5e988SBarry Smith } 114ede5e988SBarry Smith 115ede5e988SBarry Smith /* form inverse(u)*inverse(l) */ 11625783f72SBarry Smith for (kb = 1; kb <= 2; ++kb) { 11725783f72SBarry Smith k = 3 - kb; 11873d4a2d6SBarry Smith k3 = 3*k; 119ede5e988SBarry Smith kp1 = k + 1; 12073d4a2d6SBarry Smith aa = a + k3; 12125783f72SBarry Smith for (i = kp1; i <= 3; ++i) { 122b48ee343SBarry Smith work[i-1] = aa[i]; 12373d4a2d6SBarry Smith aa[i] = 0.0; 124ede5e988SBarry Smith } 12525783f72SBarry Smith for (j = kp1; j <= 3; ++j) { 126b48ee343SBarry Smith stmp = work[j-1]; 12725783f72SBarry Smith ax = &a[3*j + 1]; 12873d4a2d6SBarry Smith ay = &a[k3 + 1]; 12939d66777SBarry Smith ay[0] += stmp*ax[0]; 13039d66777SBarry Smith ay[1] += stmp*ax[1]; 13139d66777SBarry Smith ay[2] += stmp*ax[2]; 132ede5e988SBarry Smith } 133da10e913SBarry Smith l = ipvt[k-1]; 134ede5e988SBarry Smith if (l != k) { 13573d4a2d6SBarry Smith ax = &a[k3 + 1]; 13625783f72SBarry Smith ay = &a[3*l + 1]; 13773d4a2d6SBarry Smith stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp; 13873d4a2d6SBarry Smith stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp; 13973d4a2d6SBarry Smith stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp; 140ede5e988SBarry Smith } 141ede5e988SBarry Smith } 1423a40ed3dSBarry Smith PetscFunctionReturn(0); 143ede5e988SBarry Smith } 144