Actual source code: dgefa3.c
1: /*
2: Inverts 3 by 3 matrix using gaussian elimination with partial pivoting.
4: Used by the sparse factorization routines in
5: src/mat/impls/baij/seq
7: This is a combination of the Linpack routines
8: dgefa() and dgedi() specialized for a size of 3.
10: */
11: #include <petscsys.h>
13: PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_3(MatScalar *a, PetscReal shift, PetscBool allowzeropivot, PetscBool *zeropivotdetected)
14: {
15: PetscInt i__2, i__3, kp1, j, k, l, ll, i, ipvt[3], kb, k3;
16: PetscInt k4, j3;
17: MatScalar *aa, *ax, *ay, work[9], stmp;
18: MatReal tmp, max;
20: PetscFunctionBegin;
21: if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
22: shift = .333 * shift * (1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[4]) + PetscAbsScalar(a[8]));
24: /* Parameter adjustments */
25: a -= 4;
27: for (k = 1; k <= 2; ++k) {
28: kp1 = k + 1;
29: k3 = 3 * k;
30: k4 = k3 + k;
32: /* find l = pivot index */
33: i__2 = 4 - 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) {
40: max = tmp;
41: l = ll + 1;
42: }
43: }
44: l += k - 1;
45: ipvt[k - 1] = l;
47: if (a[l + k3] == 0.0) {
48: if (shift == 0.0) {
49: if (allowzeropivot) {
50: PetscCall(PetscInfo(NULL, "Zero pivot, row %" PetscInt_FMT "\n", k - 1));
51: if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
52: } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row %" PetscInt_FMT, k - 1);
53: } else {
54: /* Shift is applied to single diagonal entry */
55: a[l + k3] = shift;
56: }
57: }
59: /* interchange if necessary */
60: if (l != k) {
61: stmp = a[l + k3];
62: a[l + k3] = a[k4];
63: a[k4] = stmp;
64: }
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;
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: }
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 (PetscLikely(allowzeropivot)) {
90: PetscCall(PetscInfo(NULL, "Zero pivot, row 2\n"));
91: if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
92: } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row 2");
93: }
95: /* Now form the inverse */
96: /* compute inverse(u) */
97: for (k = 1; k <= 3; ++k) {
98: k3 = 3 * k;
99: k4 = k3 + k;
100: a[k4] = 1.0 / a[k4];
101: stmp = -a[k4];
102: i__2 = k - 1;
103: aa = &a[k3 + 1];
104: for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp;
105: kp1 = k + 1;
106: if (3 < kp1) continue;
107: ax = aa;
108: for (j = kp1; j <= 3; ++j) {
109: j3 = 3 * j;
110: stmp = a[k + j3];
111: a[k + j3] = 0.0;
112: ay = &a[j3 + 1];
113: for (ll = 0; ll < k; ll++) ay[ll] += stmp * ax[ll];
114: }
115: }
117: /* form inverse(u)*inverse(l) */
118: for (kb = 1; kb <= 2; ++kb) {
119: k = 3 - kb;
120: k3 = 3 * k;
121: kp1 = k + 1;
122: aa = a + k3;
123: for (i = kp1; i <= 3; ++i) {
124: work[i - 1] = aa[i];
125: aa[i] = 0.0;
126: }
127: for (j = kp1; j <= 3; ++j) {
128: stmp = work[j - 1];
129: ax = &a[3 * j + 1];
130: ay = &a[k3 + 1];
131: ay[0] += stmp * ax[0];
132: ay[1] += stmp * ax[1];
133: ay[2] += stmp * ax[2];
134: }
135: l = ipvt[k - 1];
136: if (l != k) {
137: ax = &a[k3 + 1];
138: ay = &a[3 * l + 1];
139: stmp = ax[0];
140: ax[0] = ay[0];
141: ay[0] = stmp;
142: stmp = ax[1];
143: ax[1] = ay[1];
144: ay[1] = stmp;
145: stmp = ax[2];
146: ax[2] = ay[2];
147: ay[2] = stmp;
148: }
149: }
150: PetscFunctionReturn(PETSC_SUCCESS);
151: }