184643e36SBarry Smith #ifdef PETSC_RCS_HEADER 2*2515f8d2SSatish Balay static char vcid[] = "$Id: dgefa2.c,v 1.1 1999/03/16 22:23:56 bsmith Exp balay $"; 384643e36SBarry Smith #endif 484643e36SBarry Smith /* 584643e36SBarry Smith Inverts 2 by 2 matrix using partial pivoting. 684643e36SBarry Smith 784643e36SBarry Smith Used by the sparse factorization routines in 884643e36SBarry Smith src/mat/impls/baij/seq and src/mat/impls/bdiag/seq 984643e36SBarry Smith 1084643e36SBarry Smith See also src/inline/ilu.h 1184643e36SBarry Smith 1284643e36SBarry Smith This is a combination of the Linpack routines 1384643e36SBarry Smith dgefa() and dgedi() specialized for a size of 2. 1484643e36SBarry Smith 1584643e36SBarry Smith */ 1684643e36SBarry Smith #include "petsc.h" 1784643e36SBarry Smith 1884643e36SBarry Smith #undef __FUNC__ 1984643e36SBarry Smith #define __FUNC__ "Kernel_A_gets_inverse_A_2" 2084643e36SBarry Smith int Kernel_A_gets_inverse_A_2(MatScalar *a) 2184643e36SBarry Smith { 22*2515f8d2SSatish Balay int i__2, i__3, kp1, j, k, l,ll,i,ipvt_l[2],*ipvt = ipvt_l-1,k3; 2384643e36SBarry Smith int k4,j3; 2484643e36SBarry Smith MatScalar *aa,*ax,*ay,work_l[4],*work = work_l-1,stmp; 2584643e36SBarry Smith MatFloat tmp,max; 2684643e36SBarry Smith 2784643e36SBarry Smith /* gaussian elimination with partial pivoting */ 2884643e36SBarry Smith 2984643e36SBarry Smith PetscFunctionBegin; 3084643e36SBarry Smith /* Parameter adjustments */ 3184643e36SBarry Smith a -= 3; 3284643e36SBarry Smith 3384643e36SBarry Smith /*for (k = 1; k <= 1; ++k) {*/ 3484643e36SBarry Smith k = 1; 3584643e36SBarry Smith kp1 = k + 1; 3684643e36SBarry Smith k3 = 2*k; 3784643e36SBarry Smith k4 = k3 + k; 3884643e36SBarry Smith /* find l = pivot index */ 3984643e36SBarry Smith 4084643e36SBarry Smith i__2 = 2 - k; 4184643e36SBarry Smith aa = &a[k4]; 4284643e36SBarry Smith max = PetscAbsScalar(aa[0]); 4384643e36SBarry Smith l = 1; 4484643e36SBarry Smith for ( ll=1; ll<i__2; ll++ ) { 4584643e36SBarry Smith tmp = PetscAbsScalar(aa[ll]); 4684643e36SBarry Smith if (tmp > max) { max = tmp; l = ll+1;} 4784643e36SBarry Smith } 4884643e36SBarry Smith l += k - 1; 4984643e36SBarry Smith ipvt[k] = l; 5084643e36SBarry Smith 5184643e36SBarry Smith if (a[l + k3] == 0.) { 5284643e36SBarry Smith SETERRQ(k,0,"Zero pivot"); 5384643e36SBarry Smith } 5484643e36SBarry Smith 5584643e36SBarry Smith /* interchange if necessary */ 5684643e36SBarry Smith 5784643e36SBarry Smith if (l != k) { 5884643e36SBarry Smith stmp = a[l + k3]; 5984643e36SBarry Smith a[l + k3] = a[k4]; 6084643e36SBarry Smith a[k4] = stmp; 6184643e36SBarry Smith } 6284643e36SBarry Smith 6384643e36SBarry Smith /* compute multipliers */ 6484643e36SBarry Smith 6584643e36SBarry Smith stmp = -1. / a[k4]; 6684643e36SBarry Smith i__2 = 2 - k; 6784643e36SBarry Smith aa = &a[1 + k4]; 6884643e36SBarry Smith for ( ll=0; ll<i__2; ll++ ) { 6984643e36SBarry Smith aa[ll] *= stmp; 7084643e36SBarry Smith } 7184643e36SBarry Smith 7284643e36SBarry Smith /* row elimination with column indexing */ 7384643e36SBarry Smith 7484643e36SBarry Smith ax = &a[k4+1]; 7584643e36SBarry Smith for (j = kp1; j <= 2; ++j) { 7684643e36SBarry Smith j3 = 2*j; 7784643e36SBarry Smith stmp = a[l + j3]; 7884643e36SBarry Smith if (l != k) { 7984643e36SBarry Smith a[l + j3] = a[k + j3]; 8084643e36SBarry Smith a[k + j3] = stmp; 8184643e36SBarry Smith } 8284643e36SBarry Smith 8384643e36SBarry Smith i__3 = 2 - k; 8484643e36SBarry Smith ay = &a[1+k+j3]; 8584643e36SBarry Smith for ( ll=0; ll<i__3; ll++ ) { 8684643e36SBarry Smith ay[ll] += stmp*ax[ll]; 8784643e36SBarry Smith } 8884643e36SBarry Smith } 8984643e36SBarry Smith /*}*/ 9084643e36SBarry Smith ipvt[2] = 2; 9184643e36SBarry Smith if (a[6] == 0.) { 9284643e36SBarry Smith SETERRQ(3,0,"Zero pivot,final row"); 9384643e36SBarry Smith } 9484643e36SBarry Smith 9584643e36SBarry Smith /* 9684643e36SBarry Smith Now form the inverse 9784643e36SBarry Smith */ 9884643e36SBarry Smith 9984643e36SBarry Smith /* compute inverse(u) */ 10084643e36SBarry Smith 10184643e36SBarry Smith for (k = 1; k <= 2; ++k) { 10284643e36SBarry Smith k3 = 2*k; 10384643e36SBarry Smith k4 = k3 + k; 10484643e36SBarry Smith a[k4] = 1.0 / a[k4]; 10584643e36SBarry Smith stmp = -a[k4]; 10684643e36SBarry Smith i__2 = k - 1; 10784643e36SBarry Smith aa = &a[k3 + 1]; 10884643e36SBarry Smith for ( ll=0; ll<i__2; ll++ ) aa[ll] *= stmp; 10984643e36SBarry Smith kp1 = k + 1; 11084643e36SBarry Smith if (2 < kp1) continue; 11184643e36SBarry Smith ax = aa; 11284643e36SBarry Smith for (j = kp1; j <= 2; ++j) { 11384643e36SBarry Smith j3 = 2*j; 11484643e36SBarry Smith stmp = a[k + j3]; 11584643e36SBarry Smith a[k + j3] = 0.0; 11684643e36SBarry Smith ay = &a[j3 + 1]; 11784643e36SBarry Smith for ( ll=0; ll<k; ll++ ) { 11884643e36SBarry Smith ay[ll] += stmp*ax[ll]; 11984643e36SBarry Smith } 12084643e36SBarry Smith } 12184643e36SBarry Smith } 12284643e36SBarry Smith 12384643e36SBarry Smith /* form inverse(u)*inverse(l) */ 12484643e36SBarry Smith 12584643e36SBarry Smith /*for (kb = 1; kb <= 1; ++kb) {*/ 12684643e36SBarry Smith 12784643e36SBarry Smith k = 1; 12884643e36SBarry Smith k3 = 2*k; 12984643e36SBarry Smith kp1 = k + 1; 13084643e36SBarry Smith aa = a + k3; 13184643e36SBarry Smith for (i = kp1; i <= 2; ++i) { 13284643e36SBarry Smith work_l[i-1] = aa[i]; 13384643e36SBarry Smith /* work[i] = aa[i]; Fix for -O3 error on Origin 2000 */ 13484643e36SBarry Smith aa[i] = 0.0; 13584643e36SBarry Smith } 13684643e36SBarry Smith for (j = kp1; j <= 2; ++j) { 13784643e36SBarry Smith stmp = work[j]; 13884643e36SBarry Smith ax = &a[2*j + 1]; 13984643e36SBarry Smith ay = &a[k3 + 1]; 14084643e36SBarry Smith ay[0] += stmp*ax[0]; 14184643e36SBarry Smith ay[1] += stmp*ax[1]; 14284643e36SBarry Smith } 14384643e36SBarry Smith l = ipvt[k]; 14484643e36SBarry Smith if (l != k) { 14584643e36SBarry Smith ax = &a[k3 + 1]; 14684643e36SBarry Smith ay = &a[2*l + 1]; 14784643e36SBarry Smith stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp; 14884643e36SBarry Smith stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp; 14984643e36SBarry Smith } 15084643e36SBarry Smith 15184643e36SBarry Smith PetscFunctionReturn(0); 15284643e36SBarry Smith } 15384643e36SBarry Smith 15484643e36SBarry Smith 155