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 PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_3(MatScalar *a,PetscReal shift,PetscBool allowzeropivot,PetscBool *zeropivotdetected) 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 *zeropivotdetected = PETSC_FALSE; 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) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1); 51 else { 52 /* Shift is applied to single diagonal entry */ 53 a[l + k3] = shift; 54 } 55 } 56 /* interchange if necessary */ 57 58 if (l != k) { 59 stmp = a[l + k3]; 60 a[l + k3] = a[k4]; 61 a[k4] = stmp; 62 } 63 64 /* compute multipliers */ 65 66 stmp = -1. / a[k4]; 67 i__2 = 3 - k; 68 aa = &a[1 + k4]; 69 for (ll=0; ll<i__2; ll++) aa[ll] *= stmp; 70 71 /* row elimination with column indexing */ 72 73 ax = &a[k4+1]; 74 for (j = kp1; j <= 3; ++j) { 75 j3 = 3*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 = 3 - k; 83 ay = &a[1+k+j3]; 84 for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll]; 85 } 86 } 87 ipvt[2] = 3; 88 if (a[12] == 0.0) { 89 PetscErrorCode ierr; 90 if (allowzeropivot) { 91 ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",2);CHKERRQ(ierr); 92 *zeropivotdetected = PETSC_TRUE; 93 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",2); 94 } 95 96 /* 97 Now form the inverse 98 */ 99 100 /* compute inverse(u) */ 101 102 for (k = 1; k <= 3; ++k) { 103 k3 = 3*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 (3 < kp1) continue; 112 ax = aa; 113 for (j = kp1; j <= 3; ++j) { 114 j3 = 3*j; 115 stmp = a[k + j3]; 116 a[k + j3] = 0.0; 117 ay = &a[j3 + 1]; 118 for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll]; 119 } 120 } 121 122 /* form inverse(u)*inverse(l) */ 123 124 for (kb = 1; kb <= 2; ++kb) { 125 k = 3 - kb; 126 k3 = 3*k; 127 kp1 = k + 1; 128 aa = a + k3; 129 for (i = kp1; i <= 3; ++i) { 130 work[i-1] = aa[i]; 131 aa[i] = 0.0; 132 } 133 for (j = kp1; j <= 3; ++j) { 134 stmp = work[j-1]; 135 ax = &a[3*j + 1]; 136 ay = &a[k3 + 1]; 137 ay[0] += stmp*ax[0]; 138 ay[1] += stmp*ax[1]; 139 ay[2] += stmp*ax[2]; 140 } 141 l = ipvt[k-1]; 142 if (l != k) { 143 ax = &a[k3 + 1]; 144 ay = &a[3*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 stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp; 148 } 149 } 150 PetscFunctionReturn(0); 151 } 152 153 154