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