Actual source code: dgefa3.c
petsc-3.7.3 2016-08-01
2: /*
3: Inverts 3 by 3 matrix using gaussian elimination with partial pivoting.
5: Used by the sparse factorization routines in
6: src/mat/impls/baij/seq
9: This is a combination of the Linpack routines
10: dgefa() and dgedi() specialized for a size of 3.
12: */
13: #include <petscsys.h>
17: PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_3(MatScalar *a,PetscReal shift,PetscBool allowzeropivot,PetscBool *zeropivotdetected)
18: {
19: PetscInt i__2,i__3,kp1,j,k,l,ll,i,ipvt[3],kb,k3;
20: PetscInt k4,j3;
21: MatScalar *aa,*ax,*ay,work[9],stmp;
22: MatReal tmp,max;
25: if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
26: shift = .333*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[4]) + PetscAbsScalar(a[8]));
28: /* Parameter adjustments */
29: a -= 4;
31: for (k = 1; k <= 2; ++k) {
32: kp1 = k + 1;
33: k3 = 3*k;
34: k4 = k3 + k;
36: /* find l = pivot index */
37: i__2 = 4 - 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: if (allowzeropivot) {
52: PetscInfo1(NULL,"Zero pivot, row %D\n",k-1);
53: if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
54: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
55: } else {
56: /* Shift is applied to single diagonal entry */
57: a[l + k3] = shift;
58: }
59: }
61: /* interchange if necessary */
62: if (l != k) {
63: stmp = a[l + k3];
64: a[l + k3] = a[k4];
65: a[k4] = stmp;
66: }
68: /* compute multipliers */
69: stmp = -1. / a[k4];
70: i__2 = 3 - k;
71: aa = &a[1 + k4];
72: for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
74: /* row elimination with column indexing */
75: ax = &a[k4+1];
76: for (j = kp1; j <= 3; ++j) {
77: j3 = 3*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 = 3 - k;
85: ay = &a[1+k+j3];
86: for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll];
87: }
88: }
89: ipvt[2] = 3;
90: if (a[12] == 0.0) {
91: if (allowzeropivot) {
93: PetscInfo1(NULL,"Zero pivot, row %D\n",2);
94: if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
95: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",2);
96: }
98: /* Now form the inverse */
99: /* compute inverse(u) */
100: for (k = 1; k <= 3; ++k) {
101: k3 = 3*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 (3 < kp1) continue;
110: ax = aa;
111: for (j = kp1; j <= 3; ++j) {
112: j3 = 3*j;
113: stmp = a[k + j3];
114: a[k + j3] = 0.0;
115: ay = &a[j3 + 1];
116: for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll];
117: }
118: }
120: /* form inverse(u)*inverse(l) */
121: for (kb = 1; kb <= 2; ++kb) {
122: k = 3 - kb;
123: k3 = 3*k;
124: kp1 = k + 1;
125: aa = a + k3;
126: for (i = kp1; i <= 3; ++i) {
127: work[i-1] = aa[i];
128: aa[i] = 0.0;
129: }
130: for (j = kp1; j <= 3; ++j) {
131: stmp = work[j-1];
132: ax = &a[3*j + 1];
133: ay = &a[k3 + 1];
134: ay[0] += stmp*ax[0];
135: ay[1] += stmp*ax[1];
136: ay[2] += stmp*ax[2];
137: }
138: l = ipvt[k-1];
139: if (l != k) {
140: ax = &a[k3 + 1];
141: ay = &a[3*l + 1];
142: stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
143: stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
144: stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
145: }
146: }
147: return(0);
148: }