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