xref: /petsc/src/mat/impls/baij/seq/dgefa7.c (revision 2e92ee13a8395f820cc1e3fd74a7607ed52efa2a)
1 
2 /*
3       Inverts 7 by 7 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 7.
10 
11 */
12 #include <petscsys.h>
13 
14 #undef __FUNCT__
15 #define __FUNCT__ "PetscKernel_A_gets_inverse_A_7"
16 PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_7(MatScalar *a,PetscReal shift,PetscBool allowzeropivot,PetscBool *zeropivotdetected)
17 {
18   PetscInt  i__2,i__3,kp1,j,k,l,ll,i,ipvt[7],kb,k3;
19   PetscInt  k4,j3;
20   MatScalar *aa,*ax,*ay,work[49],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[8]) + PetscAbsScalar(a[16]) + PetscAbsScalar(a[24]) + PetscAbsScalar(a[32]) + PetscAbsScalar(a[40]) + PetscAbsScalar(a[48]));
27 
28   /* Parameter adjustments */
29   a -= 8;
30 
31   for (k = 1; k <= 6; ++k) {
32     kp1 = k + 1;
33     k3  = 7*k;
34     k4  = k3 + k;
35     /* find l = pivot index */
36 
37     i__2 = 8 - k;
38     aa   = &a[k4];
39     max  = PetscAbsScalar(aa[0]);
40     l    = 1;
41     for (ll=1; ll<i__2; ll++) {
42       tmp = PetscAbsScalar(aa[ll]);
43       if (tmp > max) { max = tmp; l = ll+1;}
44     }
45     l        += k - 1;
46     ipvt[k-1] = l;
47 
48     if (a[l + k3] == 0.0) {
49       if (shift == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
50       else {
51         /* SHIFT is applied to SINGLE diagonal entry; does this make any sense? */
52         a[l + k3] = shift;
53       }
54     }
55 
56     /* interchange if necessary */
57 
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 
66     stmp = -1. / a[k4];
67     i__2 = 7 - k;
68     aa   = &a[1 + k4];
69     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
70 
71     /* row elimination with column indexing */
72 
73     ax = &a[k4+1];
74     for (j = kp1; j <= 7; ++j) {
75       j3   = 7*j;
76       stmp = a[l + j3];
77       if (l != k) {
78         a[l + j3] = a[k + j3];
79         a[k + j3] = stmp;
80       }
81 
82       i__3 = 7 - k;
83       ay   = &a[1+k+j3];
84       for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll];
85     }
86   }
87   ipvt[6] = 7;
88   if (a[56] == 0.0) {
89     PetscErrorCode ierr;
90     if (allowzeropivot) {
91       ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",6);CHKERRQ(ierr);
92       *zeropivotdetected = PETSC_TRUE;
93     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",6);
94   }
95 
96   /*
97      Now form the inverse
98   */
99 
100   /* compute inverse(u) */
101 
102   for (k = 1; k <= 7; ++k) {
103     k3    = 7*k;
104     k4    = k3 + k;
105     a[k4] = 1.0 / a[k4];
106     stmp  = -a[k4];
107     i__2  = k - 1;
108     aa    = &a[k3 + 1];
109     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
110     kp1 = k + 1;
111     if (7 < kp1) continue;
112     ax = aa;
113     for (j = kp1; j <= 7; ++j) {
114       j3        = 7*j;
115       stmp      = a[k + j3];
116       a[k + j3] = 0.0;
117       ay        = &a[j3 + 1];
118       for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll];
119     }
120   }
121 
122   /* form inverse(u)*inverse(l) */
123 
124   for (kb = 1; kb <= 6; ++kb) {
125     k   = 7 - kb;
126     k3  = 7*k;
127     kp1 = k + 1;
128     aa  = a + k3;
129     for (i = kp1; i <= 7; ++i) {
130       work[i-1] = aa[i];
131       aa[i]     = 0.0;
132     }
133     for (j = kp1; j <= 7; ++j) {
134       stmp   = work[j-1];
135       ax     = &a[7*j + 1];
136       ay     = &a[k3 + 1];
137       ay[0] += stmp*ax[0];
138       ay[1] += stmp*ax[1];
139       ay[2] += stmp*ax[2];
140       ay[3] += stmp*ax[3];
141       ay[4] += stmp*ax[4];
142       ay[5] += stmp*ax[5];
143       ay[6] += stmp*ax[6];
144     }
145     l = ipvt[k-1];
146     if (l != k) {
147       ax   = &a[k3 + 1];
148       ay   = &a[7*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       stmp = ax[6]; ax[6] = ay[6]; ay[6] = stmp;
156     }
157   }
158   PetscFunctionReturn(0);
159 }
160 
161 
162