xref: /petsc/src/mat/impls/baij/seq/dgefa6.c (revision c80103dabbafe9c3104be5e968fb1ec50fcbf95c)
1 
2 /*
3       Inverts 6 by 6 matrix using 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 6.
10 
11 */
12 #include <petscsys.h>
13 
14 #undef __FUNCT__
15 #define __FUNCT__ "PetscKernel_A_gets_inverse_A_6"
16 PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_6(MatScalar *a,PetscReal shift,PetscBool allowzeropivot,PetscBool *zeropivotdetected)
17 {
18   PetscInt  i__2,i__3,kp1,j,k,l,ll,i,ipvt[6],kb,k3;
19   PetscInt  k4,j3;
20   MatScalar *aa,*ax,*ay,work[36],stmp;
21   MatReal   tmp,max;
22 
23   /* gaussian elimination with partial pivoting */
24 
25   PetscFunctionBegin;
26   if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
27 
28   shift = .25*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[7]) + PetscAbsScalar(a[14]) + PetscAbsScalar(a[21]) + PetscAbsScalar(a[28]) + PetscAbsScalar(a[35]));
29   /* Parameter adjustments */
30   a -= 7;
31 
32   for (k = 1; k <= 5; ++k) {
33     kp1 = k + 1;
34     k3  = 6*k;
35     k4  = k3 + k;
36     /* find l = pivot index */
37 
38     i__2 = 7 - k;
39     aa   = &a[k4];
40     max  = PetscAbsScalar(aa[0]);
41     l    = 1;
42     for (ll=1; ll<i__2; ll++) {
43       tmp = PetscAbsScalar(aa[ll]);
44       if (tmp > max) { max = tmp; l = ll+1;}
45     }
46     l        += k - 1;
47     ipvt[k-1] = l;
48 
49     if (a[l + k3] == 0.0) {
50       if (shift == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
51       else {
52         /* SHIFT is applied to SINGLE diagonal entry; does this make any sense? */
53         a[l + k3] = shift;
54       }
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 = 6 - 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 <= 6; ++j) {
76       j3   = 6*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 = 6 - k;
84       ay   = &a[1+k+j3];
85       for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll];
86     }
87   }
88   ipvt[5] = 6;
89   if (a[42] == 0.0) {
90     PetscErrorCode ierr;
91     if (allowzeropivot) {
92       ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",5);CHKERRQ(ierr);
93       *zeropivotdetected = PETSC_TRUE;
94     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",5);
95   }
96 
97   /*
98    Now form the inverse
99   */
100 
101   /* compute inverse(u) */
102 
103   for (k = 1; k <= 6; ++k) {
104     k3    = 6*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 (6 < kp1) continue;
113     ax = aa;
114     for (j = kp1; j <= 6; ++j) {
115       j3        = 6*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 <= 5; ++kb) {
126     k   = 6 - kb;
127     k3  = 6*k;
128     kp1 = k + 1;
129     aa  = a + k3;
130     for (i = kp1; i <= 6; ++i) {
131       work[i-1] = aa[i];
132       aa[i]     = 0.0;
133     }
134     for (j = kp1; j <= 6; ++j) {
135       stmp   = work[j-1];
136       ax     = &a[6*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       ay[3] += stmp*ax[3];
142       ay[4] += stmp*ax[4];
143       ay[5] += stmp*ax[5];
144     }
145     l = ipvt[k-1];
146     if (l != k) {
147       ax   = &a[k3 + 1];
148       ay   = &a[6*l + 1];
149       stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
150       stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
151       stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
152       stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp;
153       stmp = ax[4]; ax[4] = ay[4]; ay[4] = stmp;
154       stmp = ax[5]; ax[5] = ay[5]; ay[5] = stmp;
155     }
156   }
157   PetscFunctionReturn(0);
158 }
159 
160