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) { 48ec1892c8SHong Zhang PetscErrorCode ierr; 49*c0aa6a63SJacob Faibussowitsch ierr = PetscInfo1(NULL,"Zero pivot, row %" PetscInt_FMT "\n",k-1);CHKERRQ(ierr); 5070d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE; 51*c0aa6a63SJacob Faibussowitsch } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %" PetscInt_FMT,k-1); 52ec1892c8SHong Zhang } else { 53ab055debSShri Abhyankar /* Shift is applied to single diagonal entry */ 54ab055debSShri Abhyankar a[l + k3] = shift; 55ab055debSShri Abhyankar } 56ab055debSShri Abhyankar } 57ede5e988SBarry Smith 58ec1892c8SHong Zhang /* interchange if necessary */ 59ede5e988SBarry Smith if (l != k) { 6039d66777SBarry Smith stmp = a[l + k3]; 6173d4a2d6SBarry Smith a[l + k3] = a[k4]; 6239d66777SBarry Smith a[k4] = stmp; 63ede5e988SBarry Smith } 64ede5e988SBarry Smith 65ede5e988SBarry Smith /* compute multipliers */ 6639d66777SBarry Smith stmp = -1. / a[k4]; 67ede5e988SBarry Smith i__2 = 3 - k; 6873d4a2d6SBarry Smith aa = &a[1 + k4]; 6926fbe8dcSKarl Rupp for (ll=0; ll<i__2; ll++) aa[ll] *= stmp; 70ede5e988SBarry Smith 71ede5e988SBarry Smith /* row elimination with column indexing */ 7273d4a2d6SBarry Smith ax = &a[k4+1]; 73ede5e988SBarry Smith for (j = kp1; j <= 3; ++j) { 7473d4a2d6SBarry Smith j3 = 3*j; 7539d66777SBarry Smith stmp = a[l + j3]; 76ede5e988SBarry Smith if (l != k) { 7773d4a2d6SBarry Smith a[l + j3] = a[k + j3]; 7839d66777SBarry Smith a[k + j3] = stmp; 79ede5e988SBarry Smith } 80ede5e988SBarry Smith 81ede5e988SBarry Smith i__3 = 3 - k; 8239d66777SBarry Smith ay = &a[1+k+j3]; 8326fbe8dcSKarl Rupp for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll]; 84ede5e988SBarry Smith } 85ede5e988SBarry Smith } 86da10e913SBarry Smith ipvt[2] = 3; 87cd545bacSHong Zhang if (a[12] == 0.0) { 88*c0aa6a63SJacob Faibussowitsch if (PetscLikely(allowzeropivot)) { 89ec1892c8SHong Zhang PetscErrorCode ierr; 90*c0aa6a63SJacob Faibussowitsch ierr = PetscInfo(NULL,"Zero pivot, row 2\n");CHKERRQ(ierr); 9170d8d27cSBarry Smith if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE; 92*c0aa6a63SJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row 2"); 93cd545bacSHong Zhang } 94ede5e988SBarry Smith 95ec1892c8SHong Zhang /* Now form the inverse */ 96ede5e988SBarry Smith /* compute inverse(u) */ 9725783f72SBarry Smith for (k = 1; k <= 3; ++k) { 9873d4a2d6SBarry Smith k3 = 3*k; 9973d4a2d6SBarry Smith k4 = k3 + k; 10073d4a2d6SBarry Smith a[k4] = 1.0 / a[k4]; 10139d66777SBarry Smith stmp = -a[k4]; 102ede5e988SBarry Smith i__2 = k - 1; 10373d4a2d6SBarry Smith aa = &a[k3 + 1]; 10439d66777SBarry Smith for (ll=0; ll<i__2; ll++) aa[ll] *= stmp; 105ede5e988SBarry Smith kp1 = k + 1; 10625783f72SBarry Smith if (3 < kp1) continue; 107ede5e988SBarry Smith ax = aa; 10825783f72SBarry Smith for (j = kp1; j <= 3; ++j) { 10973d4a2d6SBarry Smith j3 = 3*j; 11039d66777SBarry Smith stmp = a[k + j3]; 11139d66777SBarry Smith a[k + j3] = 0.0; 11273d4a2d6SBarry Smith ay = &a[j3 + 1]; 11326fbe8dcSKarl Rupp for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll]; 114ede5e988SBarry Smith } 115ede5e988SBarry Smith } 116ede5e988SBarry Smith 117ede5e988SBarry Smith /* form inverse(u)*inverse(l) */ 11825783f72SBarry Smith for (kb = 1; kb <= 2; ++kb) { 11925783f72SBarry Smith k = 3 - kb; 12073d4a2d6SBarry Smith k3 = 3*k; 121ede5e988SBarry Smith kp1 = k + 1; 12273d4a2d6SBarry Smith aa = a + k3; 12325783f72SBarry Smith for (i = kp1; i <= 3; ++i) { 124b48ee343SBarry Smith work[i-1] = aa[i]; 12573d4a2d6SBarry Smith aa[i] = 0.0; 126ede5e988SBarry Smith } 12725783f72SBarry Smith for (j = kp1; j <= 3; ++j) { 128b48ee343SBarry Smith stmp = work[j-1]; 12925783f72SBarry Smith ax = &a[3*j + 1]; 13073d4a2d6SBarry Smith ay = &a[k3 + 1]; 13139d66777SBarry Smith ay[0] += stmp*ax[0]; 13239d66777SBarry Smith ay[1] += stmp*ax[1]; 13339d66777SBarry Smith ay[2] += stmp*ax[2]; 134ede5e988SBarry Smith } 135da10e913SBarry Smith l = ipvt[k-1]; 136ede5e988SBarry Smith if (l != k) { 13773d4a2d6SBarry Smith ax = &a[k3 + 1]; 13825783f72SBarry Smith ay = &a[3*l + 1]; 13973d4a2d6SBarry Smith stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp; 14073d4a2d6SBarry Smith stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp; 14173d4a2d6SBarry Smith stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp; 142ede5e988SBarry Smith } 143ede5e988SBarry Smith } 1443a40ed3dSBarry Smith PetscFunctionReturn(0); 145ede5e988SBarry Smith } 146ede5e988SBarry Smith 147