Actual source code: dgefa5.c
petsc-3.3-p7 2013-05-11
2: /*
3: Inverts 5 by 5 matrix using 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>
16: PetscErrorCode PetscKernel_A_gets_inverse_A_5(MatScalar *a,PetscInt *ipvt,MatScalar *work,PetscReal shift)
17: {
18: PetscInt i__2,i__3,kp1,j,k,l,ll,i,kb,k3;
19: PetscInt k4,j3;
20: MatScalar *aa,*ax,*ay,stmp;
21: MatReal tmp,max;
23: /* gaussian elimination with partial pivoting */
26: shift = .25*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[6]) + PetscAbsScalar(a[12]) + PetscAbsScalar(a[18]) + PetscAbsScalar(a[24]));
27: /* Parameter adjustments */
28: a -= 6;
30: for (k = 1; k <= 4; ++k) {
31: kp1 = k + 1;
32: k3 = 5*k;
33: k4 = k3 + k;
34: /* find l = pivot index */
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;
47: if (a[l + k3] == 0.0) {
48: if (shift == 0.0) {
49: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",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 */
58: if (l != k) {
59: stmp = a[l + k3];
60: a[l + k3] = a[k4];
61: a[k4] = stmp;
62: }
64: /* compute multipliers */
66: stmp = -1. / a[k4];
67: i__2 = 5 - k;
68: aa = &a[1 + k4];
69: for (ll=0; ll<i__2; ll++) {
70: aa[ll] *= stmp;
71: }
73: /* row elimination with column indexing */
75: ax = &a[k4+1];
76: for (j = kp1; j <= 5; ++j) {
77: j3 = 5*j;
78: stmp = a[l + j3];
79: if (l != k) {
80: a[l + j3] = a[k + j3];
81: a[k + j3] = stmp;
82: }
84: i__3 = 5 - k;
85: ay = &a[1+k+j3];
86: for (ll=0; ll<i__3; ll++) {
87: ay[ll] += stmp*ax[ll];
88: }
89: }
90: }
91: ipvt[4] = 5;
92: if (a[30] == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",4);
94: /*
95: Now form the inverse
96: */
98: /* compute inverse(u) */
100: for (k = 1; k <= 5; ++k) {
101: k3 = 5*k;
102: k4 = k3 + k;
103: a[k4] = 1.0 / a[k4];
104: stmp = -a[k4];
105: i__2 = k - 1;
106: aa = &a[k3 + 1];
107: for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
108: kp1 = k + 1;
109: if (5 < kp1) continue;
110: ax = aa;
111: for (j = kp1; j <= 5; ++j) {
112: j3 = 5*j;
113: stmp = a[k + j3];
114: a[k + j3] = 0.0;
115: ay = &a[j3 + 1];
116: for (ll=0; ll<k; ll++) {
117: ay[ll] += stmp*ax[ll];
118: }
119: }
120: }
122: /* form inverse(u)*inverse(l) */
124: for (kb = 1; kb <= 4; ++kb) {
125: k = 5 - kb;
126: k3 = 5*k;
127: kp1 = k + 1;
128: aa = a + k3;
129: for (i = kp1; i <= 5; ++i) {
130: work[i-1] = aa[i];
131: aa[i] = 0.0;
132: }
133: for (j = kp1; j <= 5; ++j) {
134: stmp = work[j-1];
135: ax = &a[5*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: ay[3] += stmp*ax[3];
141: ay[4] += stmp*ax[4];
142: }
143: l = ipvt[k-1];
144: if (l != k) {
145: ax = &a[k3 + 1];
146: ay = &a[5*l + 1];
147: stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
148: stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
149: stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
150: stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp;
151: stmp = ax[4]; ax[4] = ay[4]; ay[4] = stmp;
152: }
153: }
154: return(0);
155: }