Actual source code: dgefa5.c
2: /*
3: Inverts 5 by 5 matrix using gaussian elimination with partial pivoting.
5: Used by the sparse factorization routines in
6: src/mat/impls/baij/seq
8: This is a combination of the Linpack routines
9: dgefa() and dgedi() specialized for a size of 5.
11: */
12: #include <petscsys.h>
14: PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_5(MatScalar *a,PetscInt *ipvt,MatScalar *work,PetscReal shift,PetscBool allowzeropivot,PetscBool *zeropivotdetected)
15: {
16: PetscInt i__2,i__3,kp1,j,k,l,ll,i,kb,k3;
17: PetscInt k4,j3;
18: MatScalar *aa,*ax,*ay,stmp;
19: MatReal tmp,max;
21: if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
22: shift = .25*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[6]) + PetscAbsScalar(a[12]) + PetscAbsScalar(a[18]) + PetscAbsScalar(a[24]));
24: /* Parameter adjustments */
25: a -= 6;
27: for (k = 1; k <= 4; ++k) {
28: kp1 = k + 1;
29: k3 = 5*k;
30: k4 = k3 + k;
32: /* find l = pivot index */
33: i__2 = 6 - k;
34: aa = &a[k4];
35: max = PetscAbsScalar(aa[0]);
36: l = 1;
37: for (ll=1; ll<i__2; ll++) {
38: tmp = PetscAbsScalar(aa[ll]);
39: if (tmp > max) { max = tmp; l = ll+1;}
40: }
41: l += k - 1;
42: ipvt[k-1] = l;
44: if (a[l + k3] == 0.0) {
45: if (shift == 0.0) {
46: if (allowzeropivot) {
47: PetscInfo(NULL,"Zero pivot, row %" PetscInt_FMT "\n",k-1);
48: *zeropivotdetected = PETSC_TRUE;
49: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %" PetscInt_FMT,k-1);
50: } else {
51: /* SHIFT is applied to SINGLE diagonal entry; does this make any sense? */
52: a[l + k3] = shift;
53: }
54: }
56: /* interchange if necessary */
57: if (l != k) {
58: stmp = a[l + k3];
59: a[l + k3] = a[k4];
60: a[k4] = stmp;
61: }
63: /* compute multipliers */
64: stmp = -1. / a[k4];
65: i__2 = 5 - k;
66: aa = &a[1 + k4];
67: for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
69: /* row elimination with column indexing */
70: ax = &a[k4+1];
71: for (j = kp1; j <= 5; ++j) {
72: j3 = 5*j;
73: stmp = a[l + j3];
74: if (l != k) {
75: a[l + j3] = a[k + j3];
76: a[k + j3] = stmp;
77: }
79: i__3 = 5 - k;
80: ay = &a[1+k+j3];
81: for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll];
82: }
83: }
84: ipvt[4] = 5;
85: if (a[30] == 0.0) {
86: if (PetscLikely(allowzeropivot)) {
87: PetscInfo(NULL,"Zero pivot, row 4\n");
88: *zeropivotdetected = PETSC_TRUE;
89: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row 4");
90: }
92: /* Now form the inverse */
93: /* compute inverse(u) */
94: for (k = 1; k <= 5; ++k) {
95: k3 = 5*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 (5 < kp1) continue;
104: ax = aa;
105: for (j = kp1; j <= 5; ++j) {
106: j3 = 5*j;
107: stmp = a[k + j3];
108: a[k + j3] = 0.0;
109: ay = &a[j3 + 1];
110: for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll];
111: }
112: }
114: /* form inverse(u)*inverse(l) */
115: for (kb = 1; kb <= 4; ++kb) {
116: k = 5 - kb;
117: k3 = 5*k;
118: kp1 = k + 1;
119: aa = a + k3;
120: for (i = kp1; i <= 5; ++i) {
121: work[i-1] = aa[i];
122: aa[i] = 0.0;
123: }
124: for (j = kp1; j <= 5; ++j) {
125: stmp = work[j-1];
126: ax = &a[5*j + 1];
127: ay = &a[k3 + 1];
128: ay[0] += stmp*ax[0];
129: ay[1] += stmp*ax[1];
130: ay[2] += stmp*ax[2];
131: ay[3] += stmp*ax[3];
132: ay[4] += stmp*ax[4];
133: }
134: l = ipvt[k-1];
135: if (l != k) {
136: ax = &a[k3 + 1];
137: ay = &a[5*l + 1];
138: stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
139: stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
140: stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
141: stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp;
142: stmp = ax[4]; ax[4] = ay[4]; ay[4] = stmp;
143: }
144: }
145: return 0;
146: }