xref: /petsc/src/mat/impls/baij/seq/dgefa3.c (revision 29bbc08cd5461c366aba6645924e8ff42acd1de0)
1*29bbc08cSBarry Smith /*$Id: dgefa3.c,v 1.17 2000/04/12 04:23:32 bsmith Exp bsmith $*/
2ede5e988SBarry Smith /*
325783f72SBarry Smith      Inverts 3 by 3 matrix using partial pivoting.
471c5468dSBarry Smith 
571c5468dSBarry Smith        Used by the sparse factorization routines in
671c5468dSBarry Smith      src/mat/impls/baij/seq and src/mat/impls/bdiag/seq
771c5468dSBarry Smith 
871c5468dSBarry Smith        See also src/inline/ilu.h
971c5468dSBarry Smith 
1071c5468dSBarry Smith        This is a combination of the Linpack routines
1171c5468dSBarry Smith     dgefa() and dgedi() specialized for a size of 3.
1271c5468dSBarry Smith 
13ede5e988SBarry Smith */
14ede5e988SBarry Smith #include "petsc.h"
15ede5e988SBarry Smith 
165615d1e5SSatish Balay #undef __FUNC__
17b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"Kernel_A_gets_inverse_A_3"
183f1db9ecSBarry Smith int Kernel_A_gets_inverse_A_3(MatScalar *a)
19ede5e988SBarry Smith {
2039d66777SBarry Smith     int        i__2,i__3,kp1,j,k,l,ll,i,ipvt_l[3],*ipvt = ipvt_l-1,kb,k3;
2173d4a2d6SBarry Smith     int        k4,j3;
223f1db9ecSBarry Smith     MatScalar  *aa,*ax,*ay,work_l[9],*work = work_l-1,stmp;
23329f5518SBarry Smith     MatReal    tmp,max;
24ede5e988SBarry Smith 
25ede5e988SBarry Smith /*     gaussian elimination with partial pivoting */
26ede5e988SBarry Smith 
273a40ed3dSBarry Smith     PetscFunctionBegin;
28ede5e988SBarry Smith     /* Parameter adjustments */
29ede5e988SBarry Smith     a       -= 4;
30ede5e988SBarry Smith 
31ede5e988SBarry Smith     for (k = 1; k <= 2; ++k) {
32ede5e988SBarry Smith 	kp1 = k + 1;
3373d4a2d6SBarry Smith         k3  = 3*k;
3473d4a2d6SBarry Smith         k4  = k3 + k;
35ede5e988SBarry Smith /*        find l = pivot index */
36ede5e988SBarry Smith 
37ede5e988SBarry Smith 	i__2 = 4 - k;
3873d4a2d6SBarry Smith         aa = &a[k4];
39ede5e988SBarry Smith         max = PetscAbsScalar(aa[0]);
40ede5e988SBarry Smith         l = 1;
41ede5e988SBarry Smith         for (ll=1; ll<i__2; ll++) {
42ede5e988SBarry Smith           tmp = PetscAbsScalar(aa[ll]);
43ede5e988SBarry Smith           if (tmp > max) { max = tmp; l = ll+1;}
44ede5e988SBarry Smith         }
45ede5e988SBarry Smith         l       += k - 1;
46ede5e988SBarry Smith 	ipvt[k] = l;
47ede5e988SBarry Smith 
4873d4a2d6SBarry Smith 	if (a[l + k3] == 0.) {
49*29bbc08cSBarry Smith 	  SETERRQ(k,"Zero pivot");
50ede5e988SBarry Smith 	}
51ede5e988SBarry Smith 
52ede5e988SBarry Smith /*           interchange if necessary */
53ede5e988SBarry Smith 
54ede5e988SBarry Smith 	if (l != k) {
5539d66777SBarry Smith 	  stmp      = a[l + k3];
5673d4a2d6SBarry Smith 	  a[l + k3] = a[k4];
5739d66777SBarry Smith 	  a[k4]     = stmp;
58ede5e988SBarry Smith         }
59ede5e988SBarry Smith 
60ede5e988SBarry Smith /*           compute multipliers */
61ede5e988SBarry Smith 
6239d66777SBarry Smith 	stmp = -1. / a[k4];
63ede5e988SBarry Smith 	i__2 = 3 - k;
6473d4a2d6SBarry Smith         aa = &a[1 + k4];
65ede5e988SBarry Smith         for (ll=0; ll<i__2; ll++) {
6639d66777SBarry Smith           aa[ll] *= stmp;
67ede5e988SBarry Smith         }
68ede5e988SBarry Smith 
69ede5e988SBarry Smith /*           row elimination with column indexing */
70ede5e988SBarry Smith 
7173d4a2d6SBarry Smith 	ax = &a[k4+1];
72ede5e988SBarry Smith         for (j = kp1; j <= 3; ++j) {
7373d4a2d6SBarry Smith             j3   = 3*j;
7439d66777SBarry Smith 	    stmp = a[l + j3];
75ede5e988SBarry Smith 	    if (l != k) {
7673d4a2d6SBarry Smith 	      a[l + j3] = a[k + j3];
7739d66777SBarry Smith 	      a[k + j3] = stmp;
78ede5e988SBarry Smith             }
79ede5e988SBarry Smith 
80ede5e988SBarry Smith 	    i__3 = 3 - k;
8139d66777SBarry Smith             ay = &a[1+k+j3];
82ede5e988SBarry Smith             for (ll=0; ll<i__3; ll++) {
8339d66777SBarry Smith               ay[ll] += stmp*ax[ll];
84ede5e988SBarry Smith             }
85ede5e988SBarry Smith 	}
86ede5e988SBarry Smith     }
87ede5e988SBarry Smith     ipvt[3] = 3;
88ede5e988SBarry Smith     if (a[12] == 0.) {
89*29bbc08cSBarry Smith 	SETERRQ(3,"Zero pivot,final row");
90ede5e988SBarry Smith     }
91ede5e988SBarry Smith 
92ede5e988SBarry Smith     /*
9325783f72SBarry Smith          Now form the inverse
94ede5e988SBarry Smith     */
95ede5e988SBarry Smith 
96ede5e988SBarry Smith    /*     compute inverse(u) */
97ede5e988SBarry Smith 
9825783f72SBarry Smith     for (k = 1; k <= 3; ++k) {
9973d4a2d6SBarry Smith         k3    = 3*k;
10073d4a2d6SBarry Smith         k4    = k3 + k;
10173d4a2d6SBarry Smith 	a[k4] = 1.0 / a[k4];
10239d66777SBarry Smith 	stmp  = -a[k4];
103ede5e988SBarry Smith 	i__2  = k - 1;
10473d4a2d6SBarry Smith         aa    = &a[k3 + 1];
10539d66777SBarry Smith         for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
106ede5e988SBarry Smith 	kp1 = k + 1;
10725783f72SBarry Smith 	if (3 < kp1) continue;
108ede5e988SBarry Smith         ax = aa;
10925783f72SBarry Smith         for (j = kp1; j <= 3; ++j) {
11073d4a2d6SBarry Smith             j3        = 3*j;
11139d66777SBarry Smith 	    stmp      = a[k + j3];
11239d66777SBarry Smith 	    a[k + j3] = 0.0;
11373d4a2d6SBarry Smith             ay        = &a[j3 + 1];
114ede5e988SBarry Smith             for (ll=0; ll<k; ll++) {
11539d66777SBarry Smith               ay[ll] += stmp*ax[ll];
116ede5e988SBarry Smith             }
117ede5e988SBarry Smith 	}
118ede5e988SBarry Smith     }
119ede5e988SBarry Smith 
120ede5e988SBarry Smith    /*    form inverse(u)*inverse(l) */
121ede5e988SBarry Smith 
12225783f72SBarry Smith     for (kb = 1; kb <= 2; ++kb) {
12325783f72SBarry Smith 	k   = 3 - kb;
12473d4a2d6SBarry Smith         k3  = 3*k;
125ede5e988SBarry Smith 	kp1 = k + 1;
12673d4a2d6SBarry Smith         aa  = a + k3;
12725783f72SBarry Smith 	for (i = kp1; i <= 3; ++i) {
128f3f07278SSatish Balay             work_l[i-1] = aa[i];
129f3f07278SSatish Balay             /* work[i] = aa[i]; Fix for -O3 error on Origin 2000 */
13073d4a2d6SBarry Smith 	    aa[i]   = 0.0;
131ede5e988SBarry Smith 	}
13225783f72SBarry Smith 	for (j = kp1; j <= 3; ++j) {
13339d66777SBarry Smith 	    stmp  = work[j];
13425783f72SBarry Smith             ax    = &a[3*j + 1];
13573d4a2d6SBarry Smith             ay    = &a[k3 + 1];
13639d66777SBarry Smith             ay[0] += stmp*ax[0];
13739d66777SBarry Smith             ay[1] += stmp*ax[1];
13839d66777SBarry Smith             ay[2] += stmp*ax[2];
139ede5e988SBarry Smith 	}
140ede5e988SBarry Smith 	l = ipvt[k];
141ede5e988SBarry Smith 	if (l != k) {
14273d4a2d6SBarry Smith             ax = &a[k3 + 1];
14325783f72SBarry Smith             ay = &a[3*l + 1];
14473d4a2d6SBarry Smith             stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
14573d4a2d6SBarry Smith             stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
14673d4a2d6SBarry Smith             stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
147ede5e988SBarry Smith 	}
148ede5e988SBarry Smith     }
1493a40ed3dSBarry Smith     PetscFunctionReturn(0);
150ede5e988SBarry Smith }
151ede5e988SBarry Smith 
15271c5468dSBarry Smith 
153