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 /* 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) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1); 49 50 /* interchange if necessary */ 51 52 if (l != k) { 53 stmp = a[l + k3]; 54 a[l + k3] = a[k4]; 55 a[k4] = stmp; 56 } 57 58 /* compute multipliers */ 59 60 stmp = -1. / a[k4]; 61 i__2 = 3 - k; 62 aa = &a[1 + k4]; 63 for (ll=0; ll<i__2; ll++) { 64 aa[ll] *= stmp; 65 } 66 67 /* row elimination with column indexing */ 68 69 ax = &a[k4+1]; 70 for (j = kp1; j <= 3; ++j) { 71 j3 = 3*j; 72 stmp = a[l + j3]; 73 if (l != k) { 74 a[l + j3] = a[k + j3]; 75 a[k + j3] = stmp; 76 } 77 78 i__3 = 3 - k; 79 ay = &a[1+k+j3]; 80 for (ll=0; ll<i__3; ll++) { 81 ay[ll] += stmp*ax[ll]; 82 } 83 } 84 } 85 ipvt[2] = 3; 86 if (a[12] == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",2); 87 88 /* 89 Now form the inverse 90 */ 91 92 /* compute inverse(u) */ 93 94 for (k = 1; k <= 3; ++k) { 95 k3 = 3*k; 96 k4 = k3 + k; 97 a[k4] = 1.0 / a[k4]; 98 stmp = -a[k4]; 99 i__2 = k - 1; 100 aa = &a[k3 + 1]; 101 for (ll=0; ll<i__2; ll++) aa[ll] *= stmp; 102 kp1 = k + 1; 103 if (3 < kp1) continue; 104 ax = aa; 105 for (j = kp1; j <= 3; ++j) { 106 j3 = 3*j; 107 stmp = a[k + j3]; 108 a[k + j3] = 0.0; 109 ay = &a[j3 + 1]; 110 for (ll=0; ll<k; ll++) { 111 ay[ll] += stmp*ax[ll]; 112 } 113 } 114 } 115 116 /* form inverse(u)*inverse(l) */ 117 118 for (kb = 1; kb <= 2; ++kb) { 119 k = 3 - kb; 120 k3 = 3*k; 121 kp1 = k + 1; 122 aa = a + k3; 123 for (i = kp1; i <= 3; ++i) { 124 work[i-1] = aa[i]; 125 aa[i] = 0.0; 126 } 127 for (j = kp1; j <= 3; ++j) { 128 stmp = work[j-1]; 129 ax = &a[3*j + 1]; 130 ay = &a[k3 + 1]; 131 ay[0] += stmp*ax[0]; 132 ay[1] += stmp*ax[1]; 133 ay[2] += stmp*ax[2]; 134 } 135 l = ipvt[k-1]; 136 if (l != k) { 137 ax = &a[k3 + 1]; 138 ay = &a[3*l + 1]; 139 stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp; 140 stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp; 141 stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp; 142 } 143 } 144 PetscFunctionReturn(0); 145 } 146 147 148