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