1 2 /* 3 Inverts 3 by 3 matrix using gaussian elimination with 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 PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_3(MatScalar *a,PetscReal shift,PetscBool allowzeropivot,PetscBool *zeropivotdetected) 16 { 17 PetscInt i__2,i__3,kp1,j,k,l,ll,i,ipvt[3],kb,k3; 18 PetscInt k4,j3; 19 MatScalar *aa,*ax,*ay,work[9],stmp; 20 MatReal tmp,max; 21 22 PetscFunctionBegin; 23 if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE; 24 shift = .333*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[4]) + PetscAbsScalar(a[8])); 25 26 /* Parameter adjustments */ 27 a -= 4; 28 29 for (k = 1; k <= 2; ++k) { 30 kp1 = k + 1; 31 k3 = 3*k; 32 k4 = k3 + k; 33 34 /* find l = pivot index */ 35 i__2 = 4 - k; 36 aa = &a[k4]; 37 max = PetscAbsScalar(aa[0]); 38 l = 1; 39 for (ll=1; ll<i__2; ll++) { 40 tmp = PetscAbsScalar(aa[ll]); 41 if (tmp > max) { max = tmp; l = ll+1;} 42 } 43 l += k - 1; 44 ipvt[k-1] = l; 45 46 if (a[l + k3] == 0.0) { 47 if (shift == 0.0) { 48 if (allowzeropivot) { 49 PetscErrorCode ierr; 50 ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",k-1);CHKERRQ(ierr); 51 if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE; 52 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1); 53 } else { 54 /* Shift is applied to single diagonal entry */ 55 a[l + k3] = shift; 56 } 57 } 58 59 /* interchange if necessary */ 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 stmp = -1. / a[k4]; 68 i__2 = 3 - k; 69 aa = &a[1 + k4]; 70 for (ll=0; ll<i__2; ll++) aa[ll] *= stmp; 71 72 /* row elimination with column indexing */ 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 if (allowzeropivot) { 90 PetscErrorCode ierr; 91 ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",2);CHKERRQ(ierr); 92 if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE; 93 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",2); 94 } 95 96 /* Now form the inverse */ 97 /* compute inverse(u) */ 98 for (k = 1; k <= 3; ++k) { 99 k3 = 3*k; 100 k4 = k3 + k; 101 a[k4] = 1.0 / a[k4]; 102 stmp = -a[k4]; 103 i__2 = k - 1; 104 aa = &a[k3 + 1]; 105 for (ll=0; ll<i__2; ll++) aa[ll] *= stmp; 106 kp1 = k + 1; 107 if (3 < kp1) continue; 108 ax = aa; 109 for (j = kp1; j <= 3; ++j) { 110 j3 = 3*j; 111 stmp = a[k + j3]; 112 a[k + j3] = 0.0; 113 ay = &a[j3 + 1]; 114 for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll]; 115 } 116 } 117 118 /* form inverse(u)*inverse(l) */ 119 for (kb = 1; kb <= 2; ++kb) { 120 k = 3 - kb; 121 k3 = 3*k; 122 kp1 = k + 1; 123 aa = a + k3; 124 for (i = kp1; i <= 3; ++i) { 125 work[i-1] = aa[i]; 126 aa[i] = 0.0; 127 } 128 for (j = kp1; j <= 3; ++j) { 129 stmp = work[j-1]; 130 ax = &a[3*j + 1]; 131 ay = &a[k3 + 1]; 132 ay[0] += stmp*ax[0]; 133 ay[1] += stmp*ax[1]; 134 ay[2] += stmp*ax[2]; 135 } 136 l = ipvt[k-1]; 137 if (l != k) { 138 ax = &a[k3 + 1]; 139 ay = &a[3*l + 1]; 140 stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp; 141 stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp; 142 stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp; 143 } 144 } 145 PetscFunctionReturn(0); 146 } 147 148 149