Actual source code: dgefa7.c
petsc-3.3-p7 2013-05-11
2: /*
3: Inverts 7 by 7 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 7.
11: */
12: #include <petscsys.h>
16: PetscErrorCode PetscKernel_A_gets_inverse_A_7(MatScalar *a,PetscReal shift)
17: {
18: PetscInt i__2,i__3,kp1,j,k,l,ll,i,ipvt[7],kb,k3;
19: PetscInt k4,j3;
20: MatScalar *aa,*ax,*ay,work[49],stmp;
21: MatReal tmp,max;
23: /* gaussian elimination with partial pivoting */
26: shift = .25*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[8]) + PetscAbsScalar(a[16]) + PetscAbsScalar(a[24]) + PetscAbsScalar(a[32]) + PetscAbsScalar(a[40]) + PetscAbsScalar(a[48]));
28: /* Parameter adjustments */
29: a -= 8;
31: for (k = 1; k <= 6; ++k) {
32: kp1 = k + 1;
33: k3 = 7*k;
34: k4 = k3 + k;
35: /* find l = pivot index */
37: i__2 = 8 - k;
38: aa = &a[k4];
39: max = PetscAbsScalar(aa[0]);
40: l = 1;
41: for (ll=1; ll<i__2; ll++) {
42: tmp = PetscAbsScalar(aa[ll]);
43: if (tmp > max) { max = tmp; l = ll+1;}
44: }
45: l += k - 1;
46: ipvt[k-1] = l;
48: if (a[l + k3] == 0.0) {
49: if (shift == 0.0) {
50: 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; does this make any sense? */
53: a[l + k3] = shift;
54: }
55: }
57: /* interchange if necessary */
59: if (l != k) {
60: stmp = a[l + k3];
61: a[l + k3] = a[k4];
62: a[k4] = stmp;
63: }
65: /* compute multipliers */
67: stmp = -1. / a[k4];
68: i__2 = 7 - k;
69: aa = &a[1 + k4];
70: for (ll=0; ll<i__2; ll++) {
71: aa[ll] *= stmp;
72: }
74: /* row elimination with column indexing */
76: ax = &a[k4+1];
77: for (j = kp1; j <= 7; ++j) {
78: j3 = 7*j;
79: stmp = a[l + j3];
80: if (l != k) {
81: a[l + j3] = a[k + j3];
82: a[k + j3] = stmp;
83: }
85: i__3 = 7 - k;
86: ay = &a[1+k+j3];
87: for (ll=0; ll<i__3; ll++) {
88: ay[ll] += stmp*ax[ll];
89: }
90: }
91: }
92: ipvt[6] = 7;
93: if (a[56] == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",6);
95: /*
96: Now form the inverse
97: */
99: /* compute inverse(u) */
101: for (k = 1; k <= 7; ++k) {
102: k3 = 7*k;
103: k4 = k3 + k;
104: a[k4] = 1.0 / a[k4];
105: stmp = -a[k4];
106: i__2 = k - 1;
107: aa = &a[k3 + 1];
108: for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
109: kp1 = k + 1;
110: if (7 < kp1) continue;
111: ax = aa;
112: for (j = kp1; j <= 7; ++j) {
113: j3 = 7*j;
114: stmp = a[k + j3];
115: a[k + j3] = 0.0;
116: ay = &a[j3 + 1];
117: for (ll=0; ll<k; ll++) {
118: ay[ll] += stmp*ax[ll];
119: }
120: }
121: }
123: /* form inverse(u)*inverse(l) */
125: for (kb = 1; kb <= 6; ++kb) {
126: k = 7 - kb;
127: k3 = 7*k;
128: kp1 = k + 1;
129: aa = a + k3;
130: for (i = kp1; i <= 7; ++i) {
131: work[i-1] = aa[i];
132: aa[i] = 0.0;
133: }
134: for (j = kp1; j <= 7; ++j) {
135: stmp = work[j-1];
136: ax = &a[7*j + 1];
137: ay = &a[k3 + 1];
138: ay[0] += stmp*ax[0];
139: ay[1] += stmp*ax[1];
140: ay[2] += stmp*ax[2];
141: ay[3] += stmp*ax[3];
142: ay[4] += stmp*ax[4];
143: ay[5] += stmp*ax[5];
144: ay[6] += stmp*ax[6];
145: }
146: l = ipvt[k-1];
147: if (l != k) {
148: ax = &a[k3 + 1];
149: ay = &a[7*l + 1];
150: stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
151: stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
152: stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
153: stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp;
154: stmp = ax[4]; ax[4] = ay[4]; ay[4] = stmp;
155: stmp = ax[5]; ax[5] = ay[5]; ay[5] = stmp;
156: stmp = ax[6]; ax[6] = ay[6]; ay[6] = stmp;
157: }
158: }
159: return(0);
160: }