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