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