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