xref: /petsc/src/mat/impls/baij/seq/dgefa7.c (revision c8a8475e04bcaa43590892a5c3e60c6f87bc31f7)
1 /*$Id: dgefa7.c,v 1.10 2001/04/07 15:51:54 bsmith Exp $*/
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 and src/mat/impls/bdiag/seq
7 
8        See also src/inline/ilu.h
9 
10        This is a combination of the Linpack routines
11     dgefa() and dgedi() specialized for a size of 7.
12 
13 */
14 #include "petsc.h"
15 
16 #undef __FUNCT__
17 #define __FUNCT__ "Kernel_A_gets_inverse_A_7"
18 int Kernel_A_gets_inverse_A_7(MatScalar *a)
19 {
20     int        i__2,i__3,kp1,j,k,l,ll,i,ipvt[7],kb,k3;
21     int        k4,j3;
22     MatScalar  *aa,*ax,*ay,work[49],stmp;
23     MatReal    tmp,max;
24 
25 /*     gaussian elimination with partial pivoting */
26 
27     PetscFunctionBegin;
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 = 7 - 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.) {
49 	  SETERRQ(k,"Zero pivot");
50 	}
51 
52 /*           interchange if necessary */
53 
54 	if (l != k) {
55 	  stmp      = a[l + k3];
56 	  a[l + k3] = a[k4];
57 	  a[k4]     = stmp;
58         }
59 
60 /*           compute multipliers */
61 
62 	stmp = -1. / a[k4];
63 	i__2 = 7 - k;
64         aa = &a[1 + k4];
65         for (ll=0; ll<i__2; ll++) {
66           aa[ll] *= stmp;
67         }
68 
69 /*           row elimination with column indexing */
70 
71 	ax = &a[k4+1];
72         for (j = kp1; j <= 7; ++j) {
73             j3   = 7*j;
74 	    stmp = a[l + j3];
75 	    if (l != k) {
76 	      a[l + j3] = a[k + j3];
77 	      a[k + j3] = stmp;
78             }
79 
80 	    i__3 = 7 - k;
81             ay = &a[1+k+j3];
82             for (ll=0; ll<i__3; ll++) {
83               ay[ll] += stmp*ax[ll];
84             }
85 	}
86     }
87     ipvt[6] = 7;
88     if (a[56] == 0.) {
89 	SETERRQ(3,"Zero pivot,final row");
90     }
91 
92     /*
93          Now form the inverse
94     */
95 
96    /*     compute inverse(u) */
97 
98     for (k = 1; k <= 7; ++k) {
99         k3    = 7*k;
100         k4    = k3 + k;
101 	a[k4] = 1.0 / a[k4];
102 	stmp  = -a[k4];
103 	i__2  = k - 1;
104         aa    = &a[k3 + 1];
105         for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
106 	kp1 = k + 1;
107 	if (7 < kp1) continue;
108         ax = aa;
109         for (j = kp1; j <= 7; ++j) {
110             j3        = 7*j;
111 	    stmp      = a[k + j3];
112 	    a[k + j3] = 0.0;
113             ay        = &a[j3 + 1];
114             for (ll=0; ll<k; ll++) {
115               ay[ll] += stmp*ax[ll];
116             }
117 	}
118     }
119 
120    /*    form inverse(u)*inverse(l) */
121 
122     for (kb = 1; kb <= 6; ++kb) {
123 	k   = 7 - kb;
124         k3  = 7*k;
125 	kp1 = k + 1;
126         aa  = a + k3;
127 	for (i = kp1; i <= 7; ++i) {
128             work[i-1] = aa[i];
129 	    aa[i]   = 0.0;
130 	}
131 	for (j = kp1; j <= 7; ++j) {
132 	    stmp  = work[j-1];
133             ax    = &a[7*j + 1];
134             ay    = &a[k3 + 1];
135             ay[0] += stmp*ax[0];
136             ay[1] += stmp*ax[1];
137             ay[2] += stmp*ax[2];
138             ay[3] += stmp*ax[3];
139             ay[4] += stmp*ax[4];
140             ay[5] += stmp*ax[5];
141             ay[6] += stmp*ax[6];
142 	}
143 	l = ipvt[k-1];
144 	if (l != k) {
145             ax = &a[k3 + 1];
146             ay = &a[7*l + 1];
147             stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
148             stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
149             stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
150             stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp;
151             stmp = ax[4]; ax[4] = ay[4]; ay[4] = stmp;
152             stmp = ax[5]; ax[5] = ay[5]; ay[5] = stmp;
153             stmp = ax[6]; ax[6] = ay[6]; ay[6] = stmp;
154 	}
155     }
156     PetscFunctionReturn(0);
157 }
158 
159 
160