1 /* 2 Inverts 2 by 2 matrix using partial pivoting. 3 4 Used by the sparse factorization routines in 5 src/mat/impls/baij/seq and src/mat/impls/bdiag/seq 6 7 See also src/inline/ilu.h 8 9 This is a combination of the Linpack routines 10 dgefa() and dgedi() specialized for a size of 2. 11 12 */ 13 #include "petsc.h" 14 15 #undef __FUNCT__ 16 #define __FUNCT__ "Kernel_A_gets_inverse_A_2" 17 PetscErrorCode Kernel_A_gets_inverse_A_2(MatScalar *a) 18 { 19 PetscInt i__2,i__3,kp1,j,k,l,ll,i,ipvt[2],k3; 20 PetscInt k4,j3; 21 MatScalar *aa,*ax,*ay,work[4],stmp; 22 MatReal tmp,max; 23 24 /* gaussian elimination with partial pivoting */ 25 26 PetscFunctionBegin; 27 /* Parameter adjustments */ 28 a -= 3; 29 30 /*for (k = 1; k <= 1; ++k) {*/ 31 k = 1; 32 kp1 = k + 1; 33 k3 = 2*k; 34 k4 = k3 + k; 35 /* find l = pivot index */ 36 37 i__2 = 2 - k; 38 aa = &a[k4]; 39 max = PetscAbsScalar(aa[0]); 40 l = 1; 41 for (ll=1; ll<i__2; ll++) { 42 tmp = PetscAbsScalar(aa[ll]); 43 if (tmp > max) { max = tmp; l = ll+1;} 44 } 45 l += k - 1; 46 ipvt[k-1] = l; 47 48 if (a[l + k3] == 0.0) { 49 SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1); 50 } 51 52 /* interchange if necessary */ 53 54 if (l != k) { 55 stmp = a[l + k3]; 56 a[l + k3] = a[k4]; 57 a[k4] = stmp; 58 } 59 60 /* compute multipliers */ 61 62 stmp = -1. / a[k4]; 63 i__2 = 2 - k; 64 aa = &a[1 + k4]; 65 for (ll=0; ll<i__2; ll++) { 66 aa[ll] *= stmp; 67 } 68 69 /* row elimination with column indexing */ 70 71 ax = &a[k4+1]; 72 for (j = kp1; j <= 2; ++j) { 73 j3 = 2*j; 74 stmp = a[l + j3]; 75 if (l != k) { 76 a[l + j3] = a[k + j3]; 77 a[k + j3] = stmp; 78 } 79 80 i__3 = 2 - k; 81 ay = &a[1+k+j3]; 82 for (ll=0; ll<i__3; ll++) { 83 ay[ll] += stmp*ax[ll]; 84 } 85 } 86 /*}*/ 87 ipvt[1] = 2; 88 if (a[6] == 0.0) { 89 SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",1); 90 } 91 92 /* 93 Now form the inverse 94 */ 95 96 /* compute inverse(u) */ 97 98 for (k = 1; k <= 2; ++k) { 99 k3 = 2*k; 100 k4 = k3 + k; 101 a[k4] = 1.0 / a[k4]; 102 stmp = -a[k4]; 103 i__2 = k - 1; 104 aa = &a[k3 + 1]; 105 for (ll=0; ll<i__2; ll++) aa[ll] *= stmp; 106 kp1 = k + 1; 107 if (2 < kp1) continue; 108 ax = aa; 109 for (j = kp1; j <= 2; ++j) { 110 j3 = 2*j; 111 stmp = a[k + j3]; 112 a[k + j3] = 0.0; 113 ay = &a[j3 + 1]; 114 for (ll=0; ll<k; ll++) { 115 ay[ll] += stmp*ax[ll]; 116 } 117 } 118 } 119 120 /* form inverse(u)*inverse(l) */ 121 122 /*for (kb = 1; kb <= 1; ++kb) {*/ 123 124 k = 1; 125 k3 = 2*k; 126 kp1 = k + 1; 127 aa = a + k3; 128 for (i = kp1; i <= 2; ++i) { 129 work[i-1] = aa[i]; 130 aa[i] = 0.0; 131 } 132 for (j = kp1; j <= 2; ++j) { 133 stmp = work[j-1]; 134 ax = &a[2*j + 1]; 135 ay = &a[k3 + 1]; 136 ay[0] += stmp*ax[0]; 137 ay[1] += stmp*ax[1]; 138 } 139 l = ipvt[k-1]; 140 if (l != k) { 141 ax = &a[k3 + 1]; 142 ay = &a[2*l + 1]; 143 stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp; 144 stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp; 145 } 146 147 PetscFunctionReturn(0); 148 } 149 150 151