xref: /petsc/src/mat/impls/baij/seq/dgefa3.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
1 
2 /*
3      Inverts 3 by 3 matrix using gaussian elimination with partial pivoting.
4 
5        Used by the sparse factorization routines in
6      src/mat/impls/baij/seq
7 
8        This is a combination of the Linpack routines
9     dgefa() and dgedi() specialized for a size of 3.
10 
11 */
12 #include <petscsys.h>
13 
14 PETSC_EXTERN 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) { max = tmp; l = ll+1;}
41     }
42     l        += k - 1;
43     ipvt[k-1] = l;
44 
45     if (a[l + k3] == 0.0) {
46       if (shift == 0.0) {
47         if (allowzeropivot) {
48           PetscCall(PetscInfo(NULL,"Zero pivot, row %" PetscInt_FMT "\n",k-1));
49           if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
50         } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %" PetscInt_FMT,k-1);
51       } else {
52         /* Shift is applied to single diagonal entry */
53         a[l + k3] = shift;
54       }
55     }
56 
57     /* interchange if necessary */
58     if (l != k) {
59       stmp      = a[l + k3];
60       a[l + k3] = a[k4];
61       a[k4]     = stmp;
62     }
63 
64     /* compute multipliers */
65     stmp = -1. / a[k4];
66     i__2 = 3 - k;
67     aa   = &a[1 + k4];
68     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
69 
70     /* row elimination with column indexing */
71     ax = &a[k4+1];
72     for (j = kp1; j <= 3; ++j) {
73       j3   = 3*j;
74       stmp = a[l + j3];
75       if (l != k) {
76         a[l + j3] = a[k + j3];
77         a[k + j3] = stmp;
78       }
79 
80       i__3 = 3 - k;
81       ay   = &a[1+k+j3];
82       for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll];
83     }
84   }
85   ipvt[2] = 3;
86   if (a[12] == 0.0) {
87     if (PetscLikely(allowzeropivot)) {
88       PetscCall(PetscInfo(NULL,"Zero pivot, row 2\n"));
89       if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
90     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row 2");
91   }
92 
93   /* Now form the inverse */
94   /* compute inverse(u) */
95   for (k = 1; k <= 3; ++k) {
96     k3    = 3*k;
97     k4    = k3 + k;
98     a[k4] = 1.0 / a[k4];
99     stmp  = -a[k4];
100     i__2  = k - 1;
101     aa    = &a[k3 + 1];
102     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
103     kp1 = k + 1;
104     if (3 < kp1) continue;
105     ax = aa;
106     for (j = kp1; j <= 3; ++j) {
107       j3        = 3*j;
108       stmp      = a[k + j3];
109       a[k + j3] = 0.0;
110       ay        = &a[j3 + 1];
111       for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll];
112     }
113   }
114 
115   /* form inverse(u)*inverse(l) */
116   for (kb = 1; kb <= 2; ++kb) {
117     k   = 3 - kb;
118     k3  = 3*k;
119     kp1 = k + 1;
120     aa  = a + k3;
121     for (i = kp1; i <= 3; ++i) {
122       work[i-1] = aa[i];
123       aa[i]     = 0.0;
124     }
125     for (j = kp1; j <= 3; ++j) {
126       stmp   = work[j-1];
127       ax     = &a[3*j + 1];
128       ay     = &a[k3 + 1];
129       ay[0] += stmp*ax[0];
130       ay[1] += stmp*ax[1];
131       ay[2] += stmp*ax[2];
132     }
133     l = ipvt[k-1];
134     if (l != k) {
135       ax   = &a[k3 + 1];
136       ay   = &a[3*l + 1];
137       stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
138       stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
139       stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
140     }
141   }
142   PetscFunctionReturn(0);
143 }
144