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