1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*71c5468dSBarry Smith static char vcid[] = "$Id: dgefa3.c,v 1.12 1998/12/17 22:10:39 bsmith Exp bsmith $"; 325783f72SBarry Smith #endif 4ede5e988SBarry Smith /* 525783f72SBarry Smith Inverts 3 by 3 matrix using partial pivoting. 6*71c5468dSBarry Smith 7*71c5468dSBarry Smith Used by the sparse factorization routines in 8*71c5468dSBarry Smith src/mat/impls/baij/seq and src/mat/impls/bdiag/seq 9*71c5468dSBarry Smith 10*71c5468dSBarry Smith See also src/inline/ilu.h 11*71c5468dSBarry Smith 12*71c5468dSBarry Smith This is a combination of the Linpack routines 13*71c5468dSBarry Smith dgefa() and dgedi() specialized for a size of 3. 14*71c5468dSBarry Smith 15ede5e988SBarry Smith */ 16ede5e988SBarry Smith #include "petsc.h" 17ede5e988SBarry Smith 185615d1e5SSatish Balay #undef __FUNC__ 195615d1e5SSatish Balay #define __FUNC__ "Kernel_A_gets_inverse_A_3" 203f1db9ecSBarry Smith int Kernel_A_gets_inverse_A_3(MatScalar *a) 21ede5e988SBarry Smith { 2239d66777SBarry Smith int i__2, i__3, kp1, j, k, l,ll,i,ipvt_l[3],*ipvt = ipvt_l-1,kb,k3; 2373d4a2d6SBarry Smith int k4,j3; 243f1db9ecSBarry Smith MatScalar *aa,*ax,*ay,work_l[9],*work = work_l-1,stmp; 253f1db9ecSBarry Smith MatFloat tmp,max; 26ede5e988SBarry Smith 27ede5e988SBarry Smith /* gaussian elimination with partial pivoting */ 28ede5e988SBarry Smith 293a40ed3dSBarry Smith PetscFunctionBegin; 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; 48ede5e988SBarry Smith ipvt[k] = l; 49ede5e988SBarry Smith 5073d4a2d6SBarry Smith if (a[l + k3] == 0.) { 51e3372554SBarry Smith SETERRQ(k,0,"Zero pivot"); 52ede5e988SBarry Smith } 53ede5e988SBarry Smith 54ede5e988SBarry Smith /* interchange if necessary */ 55ede5e988SBarry Smith 56ede5e988SBarry Smith if (l != k) { 5739d66777SBarry Smith stmp = a[l + k3]; 5873d4a2d6SBarry Smith a[l + k3] = a[k4]; 5939d66777SBarry Smith a[k4] = stmp; 60ede5e988SBarry Smith } 61ede5e988SBarry Smith 62ede5e988SBarry Smith /* compute multipliers */ 63ede5e988SBarry Smith 6439d66777SBarry Smith stmp = -1. / a[k4]; 65ede5e988SBarry Smith i__2 = 3 - k; 6673d4a2d6SBarry Smith aa = &a[1 + k4]; 67ede5e988SBarry Smith for ( ll=0; ll<i__2; ll++ ) { 6839d66777SBarry Smith aa[ll] *= stmp; 69ede5e988SBarry Smith } 70ede5e988SBarry Smith 71ede5e988SBarry Smith /* row elimination with column indexing */ 72ede5e988SBarry Smith 7373d4a2d6SBarry Smith ax = &a[k4+1]; 74ede5e988SBarry Smith for (j = kp1; j <= 3; ++j) { 7573d4a2d6SBarry Smith j3 = 3*j; 7639d66777SBarry Smith stmp = a[l + j3]; 77ede5e988SBarry Smith if (l != k) { 7873d4a2d6SBarry Smith a[l + j3] = a[k + j3]; 7939d66777SBarry Smith a[k + j3] = stmp; 80ede5e988SBarry Smith } 81ede5e988SBarry Smith 82ede5e988SBarry Smith i__3 = 3 - k; 8339d66777SBarry Smith ay = &a[1+k+j3]; 84ede5e988SBarry Smith for ( ll=0; ll<i__3; ll++ ) { 8539d66777SBarry Smith ay[ll] += stmp*ax[ll]; 86ede5e988SBarry Smith } 87ede5e988SBarry Smith } 88ede5e988SBarry Smith } 89ede5e988SBarry Smith ipvt[3] = 3; 90ede5e988SBarry Smith if (a[12] == 0.) { 91e3372554SBarry Smith SETERRQ(3,0,"Zero pivot,final row"); 92ede5e988SBarry Smith } 93ede5e988SBarry Smith 94ede5e988SBarry Smith /* 9525783f72SBarry Smith Now form the inverse 96ede5e988SBarry Smith */ 97ede5e988SBarry Smith 98ede5e988SBarry Smith /* compute inverse(u) */ 99ede5e988SBarry Smith 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]; 116ede5e988SBarry Smith for ( ll=0; ll<k; ll++ ) { 11739d66777SBarry Smith ay[ll] += stmp*ax[ll]; 118ede5e988SBarry Smith } 119ede5e988SBarry Smith } 120ede5e988SBarry Smith } 121ede5e988SBarry Smith 122ede5e988SBarry Smith /* form inverse(u)*inverse(l) */ 123ede5e988SBarry Smith 12425783f72SBarry Smith for (kb = 1; kb <= 2; ++kb) { 12525783f72SBarry Smith k = 3 - kb; 12673d4a2d6SBarry Smith k3 = 3*k; 127ede5e988SBarry Smith kp1 = k + 1; 12873d4a2d6SBarry Smith aa = a + k3; 12925783f72SBarry Smith for (i = kp1; i <= 3; ++i) { 130f3f07278SSatish Balay work_l[i-1] = aa[i]; 131f3f07278SSatish Balay /* work[i] = aa[i]; Fix for -O3 error on Origin 2000 */ 13273d4a2d6SBarry Smith aa[i] = 0.0; 133ede5e988SBarry Smith } 13425783f72SBarry Smith for (j = kp1; j <= 3; ++j) { 13539d66777SBarry Smith stmp = work[j]; 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 } 142ede5e988SBarry Smith l = ipvt[k]; 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 154*71c5468dSBarry Smith 155