1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 3ede5e988SBarry Smith /* 425783f72SBarry Smith Inverts 3 by 3 matrix using partial pivoting. 571c5468dSBarry Smith 671c5468dSBarry Smith Used by the sparse factorization routines in 7dd882469SBarry Smith src/mat/impls/baij/seq 871c5468dSBarry Smith 971c5468dSBarry Smith 1071c5468dSBarry Smith This is a combination of the Linpack routines 1171c5468dSBarry Smith dgefa() and dgedi() specialized for a size of 3. 1271c5468dSBarry Smith 13ede5e988SBarry Smith */ 14d382aafbSBarry Smith #include "petscsys.h" 15ede5e988SBarry Smith 164a2ae208SSatish Balay #undef __FUNCT__ 174a2ae208SSatish Balay #define __FUNCT__ "Kernel_A_gets_inverse_A_3" 1862bba022SBarry Smith PetscErrorCode Kernel_A_gets_inverse_A_3(MatScalar *a,PetscReal Shift) 19ede5e988SBarry Smith { 20690b6cddSBarry Smith PetscInt i__2,i__3,kp1,j,k,l,ll,i,ipvt[3],kb,k3; 21690b6cddSBarry Smith PetscInt k4,j3; 22b48ee343SBarry Smith MatScalar *aa,*ax,*ay,work[9],stmp; 23329f5518SBarry Smith MatReal tmp,max; 24ede5e988SBarry Smith 25ede5e988SBarry Smith /* gaussian elimination with partial pivoting */ 26ede5e988SBarry Smith 273a40ed3dSBarry Smith PetscFunctionBegin; 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 /* find l = pivot index */ 36ede5e988SBarry Smith 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 48*65e19b50SBarry Smith if (a[l + k3] == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1); 49ede5e988SBarry Smith 50ede5e988SBarry Smith /* interchange if necessary */ 51ede5e988SBarry Smith 52ede5e988SBarry Smith if (l != k) { 5339d66777SBarry Smith stmp = a[l + k3]; 5473d4a2d6SBarry Smith a[l + k3] = a[k4]; 5539d66777SBarry Smith a[k4] = stmp; 56ede5e988SBarry Smith } 57ede5e988SBarry Smith 58ede5e988SBarry Smith /* compute multipliers */ 59ede5e988SBarry Smith 6039d66777SBarry Smith stmp = -1. / a[k4]; 61ede5e988SBarry Smith i__2 = 3 - k; 6273d4a2d6SBarry Smith aa = &a[1 + k4]; 63ede5e988SBarry Smith for (ll=0; ll<i__2; ll++) { 6439d66777SBarry Smith aa[ll] *= stmp; 65ede5e988SBarry Smith } 66ede5e988SBarry Smith 67ede5e988SBarry Smith /* row elimination with column indexing */ 68ede5e988SBarry Smith 6973d4a2d6SBarry Smith ax = &a[k4+1]; 70ede5e988SBarry Smith for (j = kp1; j <= 3; ++j) { 7173d4a2d6SBarry Smith j3 = 3*j; 7239d66777SBarry Smith stmp = a[l + j3]; 73ede5e988SBarry Smith if (l != k) { 7473d4a2d6SBarry Smith a[l + j3] = a[k + j3]; 7539d66777SBarry Smith a[k + j3] = stmp; 76ede5e988SBarry Smith } 77ede5e988SBarry Smith 78ede5e988SBarry Smith i__3 = 3 - k; 7939d66777SBarry Smith ay = &a[1+k+j3]; 80ede5e988SBarry Smith for (ll=0; ll<i__3; ll++) { 8139d66777SBarry Smith ay[ll] += stmp*ax[ll]; 82ede5e988SBarry Smith } 83ede5e988SBarry Smith } 84ede5e988SBarry Smith } 85da10e913SBarry Smith ipvt[2] = 3; 86*65e19b50SBarry Smith if (a[12] == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",2); 87ede5e988SBarry Smith 88ede5e988SBarry Smith /* 8925783f72SBarry Smith Now form the inverse 90ede5e988SBarry Smith */ 91ede5e988SBarry Smith 92ede5e988SBarry Smith /* compute inverse(u) */ 93ede5e988SBarry Smith 9425783f72SBarry Smith for (k = 1; k <= 3; ++k) { 9573d4a2d6SBarry Smith k3 = 3*k; 9673d4a2d6SBarry Smith k4 = k3 + k; 9773d4a2d6SBarry Smith a[k4] = 1.0 / a[k4]; 9839d66777SBarry Smith stmp = -a[k4]; 99ede5e988SBarry Smith i__2 = k - 1; 10073d4a2d6SBarry Smith aa = &a[k3 + 1]; 10139d66777SBarry Smith for (ll=0; ll<i__2; ll++) aa[ll] *= stmp; 102ede5e988SBarry Smith kp1 = k + 1; 10325783f72SBarry Smith if (3 < kp1) continue; 104ede5e988SBarry Smith ax = aa; 10525783f72SBarry Smith for (j = kp1; j <= 3; ++j) { 10673d4a2d6SBarry Smith j3 = 3*j; 10739d66777SBarry Smith stmp = a[k + j3]; 10839d66777SBarry Smith a[k + j3] = 0.0; 10973d4a2d6SBarry Smith ay = &a[j3 + 1]; 110ede5e988SBarry Smith for (ll=0; ll<k; ll++) { 11139d66777SBarry Smith ay[ll] += stmp*ax[ll]; 112ede5e988SBarry Smith } 113ede5e988SBarry Smith } 114ede5e988SBarry Smith } 115ede5e988SBarry Smith 116ede5e988SBarry Smith /* form inverse(u)*inverse(l) */ 117ede5e988SBarry Smith 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 14771c5468dSBarry Smith 148