xref: /petsc/src/mat/impls/baij/seq/dgefa6.c (revision 047240e14af00aad1ef65e96f6fface8924f7f7e)
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 PetscErrorCode PetscKernel_A_gets_inverse_A_6(MatScalar *a,PetscReal shift)
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   shift = .25*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[7]) + PetscAbsScalar(a[14]) + PetscAbsScalar(a[21]) + PetscAbsScalar(a[28]) + PetscAbsScalar(a[35]));
27   /* Parameter adjustments */
28   a -= 7;
29 
30   for (k = 1; k <= 5; ++k) {
31     kp1 = k + 1;
32     k3  = 6*k;
33     k4  = k3 + k;
34 /*        find l = pivot index */
35 
36     i__2 = 7 - k;
37     aa   = &a[k4];
38     max  = PetscAbsScalar(aa[0]);
39     l    = 1;
40     for (ll=1; ll<i__2; ll++) {
41       tmp = PetscAbsScalar(aa[ll]);
42       if (tmp > max) { max = tmp; l = ll+1;}
43     }
44     l        += k - 1;
45     ipvt[k-1] = l;
46 
47     if (a[l + k3] == 0.0) {
48       if (shift == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
49       else {
50         /* SHIFT is applied to SINGLE diagonal entry; does this make any sense? */
51         a[l + k3] = shift;
52       }
53     }
54 
55 /*           interchange if necessary */
56 
57     if (l != k) {
58       stmp      = a[l + k3];
59       a[l + k3] = a[k4];
60       a[k4]     = stmp;
61     }
62 
63 /*           compute multipliers */
64 
65     stmp = -1. / a[k4];
66     i__2 = 6 - 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 
72     ax = &a[k4+1];
73     for (j = kp1; j <= 6; ++j) {
74       j3   = 6*j;
75       stmp = a[l + j3];
76       if (l != k) {
77         a[l + j3] = a[k + j3];
78         a[k + j3] = stmp;
79       }
80 
81       i__3 = 6 - k;
82       ay   = &a[1+k+j3];
83       for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll];
84     }
85   }
86   ipvt[5] = 6;
87   if (a[42] == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",5);
88 
89   /*
90        Now form the inverse
91   */
92 
93   /*     compute inverse(u) */
94 
95   for (k = 1; k <= 6; ++k) {
96     k3    = 6*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 (6 < kp1) continue;
105     ax = aa;
106     for (j = kp1; j <= 6; ++j) {
107       j3        = 6*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 
117   for (kb = 1; kb <= 5; ++kb) {
118     k   = 6 - kb;
119     k3  = 6*k;
120     kp1 = k + 1;
121     aa  = a + k3;
122     for (i = kp1; i <= 6; ++i) {
123       work[i-1] = aa[i];
124       aa[i]     = 0.0;
125     }
126     for (j = kp1; j <= 6; ++j) {
127       stmp   = work[j-1];
128       ax     = &a[6*j + 1];
129       ay     = &a[k3 + 1];
130       ay[0] += stmp*ax[0];
131       ay[1] += stmp*ax[1];
132       ay[2] += stmp*ax[2];
133       ay[3] += stmp*ax[3];
134       ay[4] += stmp*ax[4];
135       ay[5] += stmp*ax[5];
136     }
137     l = ipvt[k-1];
138     if (l != k) {
139       ax   = &a[k3 + 1];
140       ay   = &a[6*l + 1];
141       stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
142       stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
143       stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
144       stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp;
145       stmp = ax[4]; ax[4] = ay[4]; ay[4] = stmp;
146       stmp = ax[5]; ax[5] = ay[5]; ay[5] = stmp;
147     }
148   }
149   PetscFunctionReturn(0);
150 }
151 
152