xref: /petsc/src/mat/impls/baij/seq/dgefa7.c (revision b1ae27ea635b58bc3a5a567f99a72fa63199efe4)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: dgefa7.c,v 1.7 1998/12/24 04:10:08 bsmith Exp bsmith $";
3 #endif
4 /*
5       Inverts 7 by 7 matrix using partial pivoting.
6 
7        Used by the sparse factorization routines in
8      src/mat/impls/baij/seq and src/mat/impls/bdiag/seq
9 
10        See also src/inline/ilu.h
11 
12        This is a combination of the Linpack routines
13     dgefa() and dgedi() specialized for a size of 7.
14 
15 */
16 #include "petsc.h"
17 
18 #undef __FUNC__
19 #define __FUNC__ "Kernel_A_gets_inverse_A_7"
20 int Kernel_A_gets_inverse_A_7(MatScalar *a)
21 {
22     int        i__2, i__3, kp1, j, k, l,ll,i,ipvt_l[7],*ipvt = ipvt_l-1,kb,k3;
23     int        k4,j3;
24     MatScalar  *aa,*ax,*ay,work_l[49],*work = work_l-1,stmp;
25     MatFloat   tmp,max;
26 
27 /*     gaussian elimination with partial pivoting */
28 
29     PetscFunctionBegin;
30     /* Parameter adjustments */
31     a       -= 8;
32 
33     for (k = 1; k <= 6; ++k) {
34 	kp1 = k + 1;
35         k3  = 7*k;
36         k4  = k3 + k;
37 /*        find l = pivot index */
38 
39 	i__2 = 7 - 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] = l;
49 
50 	if (a[l + k3] == 0.) {
51 	  SETERRQ(k,0,"Zero pivot");
52 	}
53 
54 /*           interchange if necessary */
55 
56 	if (l != k) {
57 	  stmp      = a[l + k3];
58 	  a[l + k3] = a[k4];
59 	  a[k4]     = stmp;
60         }
61 
62 /*           compute multipliers */
63 
64 	stmp = -1. / a[k4];
65 	i__2 = 7 - k;
66         aa = &a[1 + k4];
67         for ( ll=0; ll<i__2; ll++ ) {
68           aa[ll] *= stmp;
69         }
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++ ) {
85               ay[ll] += stmp*ax[ll];
86             }
87 	}
88     }
89     ipvt[7] = 7;
90     if (a[56] == 0.) {
91 	SETERRQ(3,0,"Zero pivot,final row");
92     }
93 
94     /*
95          Now form the inverse
96     */
97 
98    /*     compute inverse(u) */
99 
100     for (k = 1; k <= 7; ++k) {
101         k3    = 7*k;
102         k4    = k3 + k;
103 	a[k4] = 1.0 / a[k4];
104 	stmp  = -a[k4];
105 	i__2  = k - 1;
106         aa    = &a[k3 + 1];
107         for ( ll=0; ll<i__2; ll++ ) aa[ll] *= stmp;
108 	kp1 = k + 1;
109 	if (7 < kp1) continue;
110         ax = aa;
111         for (j = kp1; j <= 7; ++j) {
112             j3        = 7*j;
113 	    stmp      = a[k + j3];
114 	    a[k + j3] = 0.0;
115             ay        = &a[j3 + 1];
116             for ( ll=0; ll<k; ll++ ) {
117               ay[ll] += stmp*ax[ll];
118             }
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_l[i-1] = aa[i];
131             /* work[i] = aa[i]; Fix for -O3 error on Origin 2000 */
132 	    aa[i]   = 0.0;
133 	}
134 	for (j = kp1; j <= 7; ++j) {
135 	    stmp  = work[j];
136             ax    = &a[7*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             ay[6] += stmp*ax[6];
145 	}
146 	l = ipvt[k];
147 	if (l != k) {
148             ax = &a[k3 + 1];
149             ay = &a[7*l + 1];
150             stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
151             stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
152             stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
153             stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp;
154             stmp = ax[4]; ax[4] = ay[4]; ay[4] = stmp;
155             stmp = ax[5]; ax[5] = ay[5]; ay[5] = stmp;
156             stmp = ax[6]; ax[6] = ay[6]; ay[6] = stmp;
157 	}
158     }
159     PetscFunctionReturn(0);
160 }
161 
162 
163