1be1d678aSKris Buschelman 2ede5e988SBarry Smith /* 325783f72SBarry Smith Inverts 3 by 3 matrix using 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> 14*cd545bacSHong Zhang //#include <petsc/private/matimpl.h> 15ede5e988SBarry Smith 164a2ae208SSatish Balay #undef __FUNCT__ 1796b95a6bSBarry Smith #define __FUNCT__ "PetscKernel_A_gets_inverse_A_3" 18*cd545bacSHong Zhang PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_3(MatScalar *a,PetscReal shift,PetscBool crash,PetscBool *wouldcrash) 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; 28*cd545bacSHong Zhang *wouldcrash = PETSC_FALSE; 290daa89e9SHong Zhang shift = .333*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[4]) + PetscAbsScalar(a[8])); 30ede5e988SBarry Smith /* Parameter adjustments */ 31ede5e988SBarry Smith a -= 4; 32ede5e988SBarry Smith 33ede5e988SBarry Smith for (k = 1; k <= 2; ++k) { 34ede5e988SBarry Smith kp1 = k + 1; 3573d4a2d6SBarry Smith k3 = 3*k; 3673d4a2d6SBarry Smith k4 = k3 + k; 37ede5e988SBarry Smith /* find l = pivot index */ 38ede5e988SBarry Smith 39ede5e988SBarry Smith i__2 = 4 - k; 4073d4a2d6SBarry Smith aa = &a[k4]; 41ede5e988SBarry Smith max = PetscAbsScalar(aa[0]); 42ede5e988SBarry Smith l = 1; 43ede5e988SBarry Smith for (ll=1; ll<i__2; ll++) { 44ede5e988SBarry Smith tmp = PetscAbsScalar(aa[ll]); 45ede5e988SBarry Smith if (tmp > max) { max = tmp; l = ll+1;} 46ede5e988SBarry Smith } 47ede5e988SBarry Smith l += k - 1; 48da10e913SBarry Smith ipvt[k-1] = l; 49ede5e988SBarry Smith 50ab055debSShri Abhyankar if (a[l + k3] == 0.0) { 5126fbe8dcSKarl Rupp if (shift == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1); 5226fbe8dcSKarl Rupp else { 53ab055debSShri Abhyankar /* Shift is applied to single diagonal entry */ 54ab055debSShri Abhyankar a[l + k3] = shift; 55ab055debSShri Abhyankar } 56ab055debSShri Abhyankar } 57ede5e988SBarry Smith /* interchange if necessary */ 58ede5e988SBarry Smith 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 */ 66ede5e988SBarry Smith 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 */ 73ede5e988SBarry Smith 7473d4a2d6SBarry Smith ax = &a[k4+1]; 75ede5e988SBarry Smith for (j = kp1; j <= 3; ++j) { 7673d4a2d6SBarry Smith j3 = 3*j; 7739d66777SBarry Smith stmp = a[l + j3]; 78ede5e988SBarry Smith if (l != k) { 7973d4a2d6SBarry Smith a[l + j3] = a[k + j3]; 8039d66777SBarry Smith a[k + j3] = stmp; 81ede5e988SBarry Smith } 82ede5e988SBarry Smith 83ede5e988SBarry Smith i__3 = 3 - k; 8439d66777SBarry Smith ay = &a[1+k+j3]; 8526fbe8dcSKarl Rupp for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll]; 86ede5e988SBarry Smith } 87ede5e988SBarry Smith } 88da10e913SBarry Smith ipvt[2] = 3; 89*cd545bacSHong Zhang if (a[12] == 0.0) { 90*cd545bacSHong Zhang PetscErrorCode ierr; 91*cd545bacSHong Zhang if (!crash) { 92*cd545bacSHong Zhang ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",2);CHKERRQ(ierr); 93*cd545bacSHong Zhang *wouldcrash = PETSC_TRUE; 94*cd545bacSHong Zhang } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",2); 95*cd545bacSHong Zhang } 96ede5e988SBarry Smith 97ede5e988SBarry Smith /* 9825783f72SBarry Smith Now form the inverse 99ede5e988SBarry Smith */ 100ede5e988SBarry Smith 101ede5e988SBarry Smith /* compute inverse(u) */ 102ede5e988SBarry Smith 10325783f72SBarry Smith for (k = 1; k <= 3; ++k) { 10473d4a2d6SBarry Smith k3 = 3*k; 10573d4a2d6SBarry Smith k4 = k3 + k; 10673d4a2d6SBarry Smith a[k4] = 1.0 / a[k4]; 10739d66777SBarry Smith stmp = -a[k4]; 108ede5e988SBarry Smith i__2 = k - 1; 10973d4a2d6SBarry Smith aa = &a[k3 + 1]; 11039d66777SBarry Smith for (ll=0; ll<i__2; ll++) aa[ll] *= stmp; 111ede5e988SBarry Smith kp1 = k + 1; 11225783f72SBarry Smith if (3 < kp1) continue; 113ede5e988SBarry Smith ax = aa; 11425783f72SBarry Smith for (j = kp1; j <= 3; ++j) { 11573d4a2d6SBarry Smith j3 = 3*j; 11639d66777SBarry Smith stmp = a[k + j3]; 11739d66777SBarry Smith a[k + j3] = 0.0; 11873d4a2d6SBarry Smith ay = &a[j3 + 1]; 11926fbe8dcSKarl Rupp for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll]; 120ede5e988SBarry Smith } 121ede5e988SBarry Smith } 122ede5e988SBarry Smith 123ede5e988SBarry Smith /* form inverse(u)*inverse(l) */ 124ede5e988SBarry Smith 12525783f72SBarry Smith for (kb = 1; kb <= 2; ++kb) { 12625783f72SBarry Smith k = 3 - kb; 12773d4a2d6SBarry Smith k3 = 3*k; 128ede5e988SBarry Smith kp1 = k + 1; 12973d4a2d6SBarry Smith aa = a + k3; 13025783f72SBarry Smith for (i = kp1; i <= 3; ++i) { 131b48ee343SBarry Smith work[i-1] = aa[i]; 13273d4a2d6SBarry Smith aa[i] = 0.0; 133ede5e988SBarry Smith } 13425783f72SBarry Smith for (j = kp1; j <= 3; ++j) { 135b48ee343SBarry Smith stmp = work[j-1]; 13625783f72SBarry Smith ax = &a[3*j + 1]; 13773d4a2d6SBarry Smith ay = &a[k3 + 1]; 13839d66777SBarry Smith ay[0] += stmp*ax[0]; 13939d66777SBarry Smith ay[1] += stmp*ax[1]; 14039d66777SBarry Smith ay[2] += stmp*ax[2]; 141ede5e988SBarry Smith } 142da10e913SBarry Smith l = ipvt[k-1]; 143ede5e988SBarry Smith if (l != k) { 14473d4a2d6SBarry Smith ax = &a[k3 + 1]; 14525783f72SBarry Smith ay = &a[3*l + 1]; 14673d4a2d6SBarry Smith stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp; 14773d4a2d6SBarry Smith stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp; 14873d4a2d6SBarry Smith stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp; 149ede5e988SBarry Smith } 150ede5e988SBarry Smith } 1513a40ed3dSBarry Smith PetscFunctionReturn(0); 152ede5e988SBarry Smith } 153ede5e988SBarry Smith 15471c5468dSBarry Smith 155