xref: /petsc/src/mat/impls/baij/seq/dgefa2.c (revision ce0a2cd1da0658c2b28aad1be2e2c8e41567bece)
1 #define PETSCMAT_DLL
2 
3 /*
4      Inverts 2 by 2 matrix using partial pivoting.
5 
6        Used by the sparse factorization routines in
7      src/mat/impls/baij/seq
8 
9        See also src/inline/ilu.h
10 
11        This is a combination of the Linpack routines
12     dgefa() and dgedi() specialized for a size of 2.
13 
14 */
15 #include "petsc.h"
16 
17 #undef __FUNCT__
18 #define __FUNCT__ "Kernel_A_gets_inverse_A_2"
19 PetscErrorCode Kernel_A_gets_inverse_A_2(MatScalar *a,PetscReal shift)
20 {
21     PetscInt   i__2,i__3,kp1,j,k,l,ll,i,ipvt[2],k3;
22     PetscInt   k4,j3;
23     MatScalar  *aa,*ax,*ay,work[4],stmp;
24     MatReal    tmp,max;
25 
26 /*     gaussian elimination with partial pivoting */
27 
28     PetscFunctionBegin;
29     /* Parameter adjustments */
30     a       -= 3;
31 
32     /*for (k = 1; k <= 1; ++k) {*/
33         k   = 1;
34 	kp1 = k + 1;
35         k3  = 2*k;
36         k4  = k3 + k;
37 /*        find l = pivot index */
38 
39 	i__2 = 3 - 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-1] = l;
49 
50 	if (a[l + k3] == 0.0) {
51 	  SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
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 = 2 - 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 <= 2; ++j) {
75             j3   = 2*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 = 2 - 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[1] = 2;
90     if (a[6] == 0.0) {
91       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",1);
92     }
93 
94     /*
95          Now form the inverse
96     */
97 
98    /*     compute inverse(u) */
99 
100     for (k = 1; k <= 2; ++k) {
101         k3    = 2*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 (2 < kp1) continue;
110         ax = aa;
111         for (j = kp1; j <= 2; ++j) {
112             j3        = 2*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 <= 1; ++kb) {*/
125 
126 	k   = 1;
127         k3  = 2*k;
128 	kp1 = k + 1;
129         aa  = a + k3;
130 	for (i = kp1; i <= 2; ++i) {
131             work[i-1] = aa[i];
132 	    aa[i]   = 0.0;
133 	}
134 	for (j = kp1; j <= 2; ++j) {
135 	    stmp  = work[j-1];
136             ax    = &a[2*j + 1];
137             ay    = &a[k3 + 1];
138             ay[0] += stmp*ax[0];
139             ay[1] += stmp*ax[1];
140 	}
141 	l = ipvt[k-1];
142 	if (l != k) {
143             ax = &a[k3 + 1];
144             ay = &a[2*l + 1];
145             stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
146             stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
147 	}
148 
149     PetscFunctionReturn(0);
150 }
151 
152 #undef __FUNCT__
153 #define __FUNCT__ "Kernel_A_gets_inverse_A_9"
154 PetscErrorCode Kernel_A_gets_inverse_A_9(MatScalar *a,PetscReal shift)
155 {
156     PetscInt   i__2,i__3,kp1,j,k,l,ll,i,ipvt[9],kb,k3;
157     PetscInt   k4,j3;
158     MatScalar  *aa,*ax,*ay,work[81],stmp;
159     MatReal    tmp,max;
160 
161 /*     gaussian elimination with partial pivoting */
162 
163     PetscFunctionBegin;
164     /* Parameter adjustments */
165     a       -= 10;
166 
167     for (k = 1; k <= 8; ++k) {
168 	kp1 = k + 1;
169         k3  = 9*k;
170         k4  = k3 + k;
171 /*        find l = pivot index */
172 
173 	i__2 = 9 - k;
174         aa = &a[k4];
175         max = PetscAbsScalar(aa[0]);
176         l = 1;
177         for (ll=1; ll<i__2; ll++) {
178           tmp = PetscAbsScalar(aa[ll]);
179           if (tmp > max) { max = tmp; l = ll+1;}
180         }
181         l       += k - 1;
182 	ipvt[k-1] = l;
183 
184 	if (a[l + k3] == 0.0) {
185 	  SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
186 	}
187 
188 /*           interchange if necessary */
189 
190 	if (l != k) {
191 	  stmp      = a[l + k3];
192 	  a[l + k3] = a[k4];
193 	  a[k4]     = stmp;
194         }
195 
196 /*           compute multipliers */
197 
198 	stmp = -1. / a[k4];
199 	i__2 = 9 - k;
200         aa = &a[1 + k4];
201         for (ll=0; ll<i__2; ll++) {
202           aa[ll] *= stmp;
203         }
204 
205 /*           row elimination with column indexing */
206 
207 	ax = &a[k4+1];
208         for (j = kp1; j <= 9; ++j) {
209             j3   = 9*j;
210 	    stmp = a[l + j3];
211 	    if (l != k) {
212 	      a[l + j3] = a[k + j3];
213 	      a[k + j3] = stmp;
214             }
215 
216 	    i__3 = 9 - k;
217             ay = &a[1+k+j3];
218             for (ll=0; ll<i__3; ll++) {
219               ay[ll] += stmp*ax[ll];
220             }
221 	}
222     }
223     ipvt[8] = 9;
224     if (a[90] == 0.0) {
225       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",6);
226     }
227 
228     /*
229          Now form the inverse
230     */
231 
232    /*     compute inverse(u) */
233 
234     for (k = 1; k <= 9; ++k) {
235         k3    = 9*k;
236         k4    = k3 + k;
237 	a[k4] = 1.0 / a[k4];
238 	stmp  = -a[k4];
239 	i__2  = k - 1;
240         aa    = &a[k3 + 1];
241         for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
242 	kp1 = k + 1;
243 	if (9 < kp1) continue;
244         ax = aa;
245         for (j = kp1; j <= 9; ++j) {
246             j3        = 9*j;
247 	    stmp      = a[k + j3];
248 	    a[k + j3] = 0.0;
249             ay        = &a[j3 + 1];
250             for (ll=0; ll<k; ll++) {
251               ay[ll] += stmp*ax[ll];
252             }
253 	}
254     }
255 
256    /*    form inverse(u)*inverse(l) */
257 
258     for (kb = 1; kb <= 8; ++kb) {
259 	k   = 9 - kb;
260         k3  = 9*k;
261 	kp1 = k + 1;
262         aa  = a + k3;
263 	for (i = kp1; i <= 9; ++i) {
264             work[i-1] = aa[i];
265 	    aa[i]   = 0.0;
266 	}
267 	for (j = kp1; j <= 9; ++j) {
268 	    stmp  = work[j-1];
269             ax    = &a[9*j + 1];
270             ay    = &a[k3 + 1];
271             ay[0] += stmp*ax[0];
272             ay[1] += stmp*ax[1];
273             ay[2] += stmp*ax[2];
274             ay[3] += stmp*ax[3];
275             ay[4] += stmp*ax[4];
276             ay[5] += stmp*ax[5];
277             ay[6] += stmp*ax[6];
278             ay[7] += stmp*ax[7];
279             ay[8] += stmp*ax[8];
280 	}
281 	l = ipvt[k-1];
282 	if (l != k) {
283             ax = &a[k3 + 1];
284             ay = &a[9*l + 1];
285             stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
286             stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
287             stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
288             stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp;
289             stmp = ax[4]; ax[4] = ay[4]; ay[4] = stmp;
290             stmp = ax[5]; ax[5] = ay[5]; ay[5] = stmp;
291             stmp = ax[6]; ax[6] = ay[6]; ay[6] = stmp;
292             stmp = ax[7]; ax[7] = ay[7]; ay[7] = stmp;
293             stmp = ax[8]; ax[8] = ay[8]; ay[8] = stmp;
294 	}
295     }
296     PetscFunctionReturn(0);
297 }
298 
299