1be1d678aSKris Buschelman 2ede5e988SBarry Smith /* 3*ec1892c8SHong 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 971c5468dSBarry Smith This is a combination of the Linpack routines 1071c5468dSBarry Smith dgefa() and dgedi() specialized for a size of 3. 1171c5468dSBarry Smith 12ede5e988SBarry Smith */ 13c6db04a5SJed Brown #include <petscsys.h> 14ede5e988SBarry Smith 154a2ae208SSatish Balay #undef __FUNCT__ 1696b95a6bSBarry Smith #define __FUNCT__ "PetscKernel_A_gets_inverse_A_3" 176baedc03SHong Zhang PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_3(MatScalar *a,PetscReal shift,PetscBool allowzeropivot,PetscBool *zeropivotdetected) 18ede5e988SBarry Smith { 19690b6cddSBarry Smith PetscInt i__2,i__3,kp1,j,k,l,ll,i,ipvt[3],kb,k3; 20690b6cddSBarry Smith PetscInt k4,j3; 21b48ee343SBarry Smith MatScalar *aa,*ax,*ay,work[9],stmp; 22329f5518SBarry Smith MatReal tmp,max; 23ede5e988SBarry Smith 243a40ed3dSBarry Smith PetscFunctionBegin; 25c80103daSHong Zhang if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE; 260daa89e9SHong Zhang shift = .333*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[4]) + PetscAbsScalar(a[8])); 27*ec1892c8SHong Zhang 28ede5e988SBarry Smith /* Parameter adjustments */ 29ede5e988SBarry Smith a -= 4; 30ede5e988SBarry Smith 31ede5e988SBarry Smith for (k = 1; k <= 2; ++k) { 32ede5e988SBarry Smith kp1 = k + 1; 3373d4a2d6SBarry Smith k3 = 3*k; 3473d4a2d6SBarry Smith k4 = k3 + k; 35ede5e988SBarry Smith 36*ec1892c8SHong Zhang /* find l = pivot index */ 37ede5e988SBarry Smith i__2 = 4 - k; 3873d4a2d6SBarry Smith aa = &a[k4]; 39ede5e988SBarry Smith max = PetscAbsScalar(aa[0]); 40ede5e988SBarry Smith l = 1; 41ede5e988SBarry Smith for (ll=1; ll<i__2; ll++) { 42ede5e988SBarry Smith tmp = PetscAbsScalar(aa[ll]); 43ede5e988SBarry Smith if (tmp > max) { max = tmp; l = ll+1;} 44ede5e988SBarry Smith } 45ede5e988SBarry Smith l += k - 1; 46da10e913SBarry Smith ipvt[k-1] = l; 47ede5e988SBarry Smith 48ab055debSShri Abhyankar if (a[l + k3] == 0.0) { 49*ec1892c8SHong Zhang if (shift == 0.0) { 50*ec1892c8SHong Zhang if (allowzeropivot) { 51*ec1892c8SHong Zhang PetscErrorCode ierr; 52*ec1892c8SHong Zhang ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",k-1);CHKERRQ(ierr); 53*ec1892c8SHong Zhang *zeropivotdetected = PETSC_TRUE; 54*ec1892c8SHong Zhang } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1); 55*ec1892c8SHong Zhang } else { 56ab055debSShri Abhyankar /* Shift is applied to single diagonal entry */ 57ab055debSShri Abhyankar a[l + k3] = shift; 58ab055debSShri Abhyankar } 59ab055debSShri Abhyankar } 60ede5e988SBarry Smith 61*ec1892c8SHong Zhang /* interchange if necessary */ 62ede5e988SBarry Smith if (l != k) { 6339d66777SBarry Smith stmp = a[l + k3]; 6473d4a2d6SBarry Smith a[l + k3] = a[k4]; 6539d66777SBarry Smith a[k4] = stmp; 66ede5e988SBarry Smith } 67ede5e988SBarry Smith 68ede5e988SBarry Smith /* compute multipliers */ 6939d66777SBarry Smith stmp = -1. / a[k4]; 70ede5e988SBarry Smith i__2 = 3 - k; 7173d4a2d6SBarry Smith aa = &a[1 + k4]; 7226fbe8dcSKarl Rupp for (ll=0; ll<i__2; ll++) aa[ll] *= stmp; 73ede5e988SBarry Smith 74ede5e988SBarry Smith /* row elimination with column indexing */ 7573d4a2d6SBarry Smith ax = &a[k4+1]; 76ede5e988SBarry Smith for (j = kp1; j <= 3; ++j) { 7773d4a2d6SBarry Smith j3 = 3*j; 7839d66777SBarry Smith stmp = a[l + j3]; 79ede5e988SBarry Smith if (l != k) { 8073d4a2d6SBarry Smith a[l + j3] = a[k + j3]; 8139d66777SBarry Smith a[k + j3] = stmp; 82ede5e988SBarry Smith } 83ede5e988SBarry Smith 84ede5e988SBarry Smith i__3 = 3 - k; 8539d66777SBarry Smith ay = &a[1+k+j3]; 8626fbe8dcSKarl Rupp for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll]; 87ede5e988SBarry Smith } 88ede5e988SBarry Smith } 89da10e913SBarry Smith ipvt[2] = 3; 90cd545bacSHong Zhang if (a[12] == 0.0) { 916baedc03SHong Zhang if (allowzeropivot) { 92*ec1892c8SHong Zhang PetscErrorCode ierr; 93cd545bacSHong Zhang ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",2);CHKERRQ(ierr); 946baedc03SHong Zhang *zeropivotdetected = PETSC_TRUE; 95cd545bacSHong Zhang } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",2); 96cd545bacSHong Zhang } 97ede5e988SBarry Smith 98*ec1892c8SHong Zhang /* Now form the inverse */ 99ede5e988SBarry Smith /* compute inverse(u) */ 10025783f72SBarry Smith for (k = 1; k <= 3; ++k) { 10173d4a2d6SBarry Smith k3 = 3*k; 10273d4a2d6SBarry Smith k4 = k3 + k; 10373d4a2d6SBarry Smith a[k4] = 1.0 / a[k4]; 10439d66777SBarry Smith stmp = -a[k4]; 105ede5e988SBarry Smith i__2 = k - 1; 10673d4a2d6SBarry Smith aa = &a[k3 + 1]; 10739d66777SBarry Smith for (ll=0; ll<i__2; ll++) aa[ll] *= stmp; 108ede5e988SBarry Smith kp1 = k + 1; 10925783f72SBarry Smith if (3 < kp1) continue; 110ede5e988SBarry Smith ax = aa; 11125783f72SBarry Smith for (j = kp1; j <= 3; ++j) { 11273d4a2d6SBarry Smith j3 = 3*j; 11339d66777SBarry Smith stmp = a[k + j3]; 11439d66777SBarry Smith a[k + j3] = 0.0; 11573d4a2d6SBarry Smith ay = &a[j3 + 1]; 11626fbe8dcSKarl Rupp for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll]; 117ede5e988SBarry Smith } 118ede5e988SBarry Smith } 119ede5e988SBarry Smith 120ede5e988SBarry Smith /* form inverse(u)*inverse(l) */ 12125783f72SBarry Smith for (kb = 1; kb <= 2; ++kb) { 12225783f72SBarry Smith k = 3 - kb; 12373d4a2d6SBarry Smith k3 = 3*k; 124ede5e988SBarry Smith kp1 = k + 1; 12573d4a2d6SBarry Smith aa = a + k3; 12625783f72SBarry Smith for (i = kp1; i <= 3; ++i) { 127b48ee343SBarry Smith work[i-1] = aa[i]; 12873d4a2d6SBarry Smith aa[i] = 0.0; 129ede5e988SBarry Smith } 13025783f72SBarry Smith for (j = kp1; j <= 3; ++j) { 131b48ee343SBarry Smith stmp = work[j-1]; 13225783f72SBarry Smith ax = &a[3*j + 1]; 13373d4a2d6SBarry Smith ay = &a[k3 + 1]; 13439d66777SBarry Smith ay[0] += stmp*ax[0]; 13539d66777SBarry Smith ay[1] += stmp*ax[1]; 13639d66777SBarry Smith ay[2] += stmp*ax[2]; 137ede5e988SBarry Smith } 138da10e913SBarry Smith l = ipvt[k-1]; 139ede5e988SBarry Smith if (l != k) { 14073d4a2d6SBarry Smith ax = &a[k3 + 1]; 14125783f72SBarry Smith ay = &a[3*l + 1]; 14273d4a2d6SBarry Smith stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp; 14373d4a2d6SBarry Smith stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp; 14473d4a2d6SBarry Smith stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp; 145ede5e988SBarry Smith } 146ede5e988SBarry Smith } 1473a40ed3dSBarry Smith PetscFunctionReturn(0); 148ede5e988SBarry Smith } 149ede5e988SBarry Smith 15071c5468dSBarry Smith 151