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