xref: /petsc/src/mat/impls/baij/seq/dgefa3.c (revision be1d678a52e6eff2808b2fa31ae986cdbf03c9fe)
1*be1d678aSKris Buschelman #define PETSCMAT_DLL
2*be1d678aSKris Buschelman 
3ede5e988SBarry Smith /*
425783f72SBarry Smith      Inverts 3 by 3 matrix using partial pivoting.
571c5468dSBarry Smith 
671c5468dSBarry Smith        Used by the sparse factorization routines in
771c5468dSBarry Smith      src/mat/impls/baij/seq and src/mat/impls/bdiag/seq
871c5468dSBarry Smith 
971c5468dSBarry Smith        See also src/inline/ilu.h
1071c5468dSBarry Smith 
1171c5468dSBarry Smith        This is a combination of the Linpack routines
1271c5468dSBarry Smith     dgefa() and dgedi() specialized for a size of 3.
1371c5468dSBarry Smith 
14ede5e988SBarry Smith */
15ede5e988SBarry Smith #include "petsc.h"
16ede5e988SBarry Smith 
174a2ae208SSatish Balay #undef __FUNCT__
184a2ae208SSatish Balay #define __FUNCT__ "Kernel_A_gets_inverse_A_3"
19dfbe8321SBarry Smith PetscErrorCode Kernel_A_gets_inverse_A_3(MatScalar *a)
20ede5e988SBarry Smith {
21690b6cddSBarry Smith     PetscInt   i__2,i__3,kp1,j,k,l,ll,i,ipvt[3],kb,k3;
22690b6cddSBarry Smith     PetscInt   k4,j3;
23b48ee343SBarry Smith     MatScalar  *aa,*ax,*ay,work[9],stmp;
24329f5518SBarry Smith     MatReal    tmp,max;
25ede5e988SBarry Smith 
26ede5e988SBarry Smith /*     gaussian elimination with partial pivoting */
27ede5e988SBarry Smith 
283a40ed3dSBarry Smith     PetscFunctionBegin;
29ede5e988SBarry Smith     /* Parameter adjustments */
30ede5e988SBarry Smith     a       -= 4;
31ede5e988SBarry Smith 
32ede5e988SBarry Smith     for (k = 1; k <= 2; ++k) {
33ede5e988SBarry Smith 	kp1 = k + 1;
3473d4a2d6SBarry Smith         k3  = 3*k;
3573d4a2d6SBarry Smith         k4  = k3 + k;
36ede5e988SBarry Smith /*        find l = pivot index */
37ede5e988SBarry Smith 
38ede5e988SBarry Smith 	i__2 = 4 - k;
3973d4a2d6SBarry Smith         aa = &a[k4];
40ede5e988SBarry Smith         max = PetscAbsScalar(aa[0]);
41ede5e988SBarry Smith         l = 1;
42ede5e988SBarry Smith         for (ll=1; ll<i__2; ll++) {
43ede5e988SBarry Smith           tmp = PetscAbsScalar(aa[ll]);
44ede5e988SBarry Smith           if (tmp > max) { max = tmp; l = ll+1;}
45ede5e988SBarry Smith         }
46ede5e988SBarry Smith         l       += k - 1;
47da10e913SBarry Smith 	ipvt[k-1] = l;
48ede5e988SBarry Smith 
495b8514ebSBarry Smith 	if (a[l + k3] == 0.0) {
5077431f27SBarry Smith 	  SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
51ede5e988SBarry Smith 	}
52ede5e988SBarry Smith 
53ede5e988SBarry Smith /*           interchange if necessary */
54ede5e988SBarry Smith 
55ede5e988SBarry Smith 	if (l != k) {
5639d66777SBarry Smith 	  stmp      = a[l + k3];
5773d4a2d6SBarry Smith 	  a[l + k3] = a[k4];
5839d66777SBarry Smith 	  a[k4]     = stmp;
59ede5e988SBarry Smith         }
60ede5e988SBarry Smith 
61ede5e988SBarry Smith /*           compute multipliers */
62ede5e988SBarry Smith 
6339d66777SBarry Smith 	stmp = -1. / a[k4];
64ede5e988SBarry Smith 	i__2 = 3 - k;
6573d4a2d6SBarry Smith         aa = &a[1 + k4];
66ede5e988SBarry Smith         for (ll=0; ll<i__2; ll++) {
6739d66777SBarry Smith           aa[ll] *= stmp;
68ede5e988SBarry Smith         }
69ede5e988SBarry Smith 
70ede5e988SBarry Smith /*           row elimination with column indexing */
71ede5e988SBarry Smith 
7273d4a2d6SBarry Smith 	ax = &a[k4+1];
73ede5e988SBarry Smith         for (j = kp1; j <= 3; ++j) {
7473d4a2d6SBarry Smith             j3   = 3*j;
7539d66777SBarry Smith 	    stmp = a[l + j3];
76ede5e988SBarry Smith 	    if (l != k) {
7773d4a2d6SBarry Smith 	      a[l + j3] = a[k + j3];
7839d66777SBarry Smith 	      a[k + j3] = stmp;
79ede5e988SBarry Smith             }
80ede5e988SBarry Smith 
81ede5e988SBarry Smith 	    i__3 = 3 - k;
8239d66777SBarry Smith             ay = &a[1+k+j3];
83ede5e988SBarry Smith             for (ll=0; ll<i__3; ll++) {
8439d66777SBarry Smith               ay[ll] += stmp*ax[ll];
85ede5e988SBarry Smith             }
86ede5e988SBarry Smith 	}
87ede5e988SBarry Smith     }
88da10e913SBarry Smith     ipvt[2] = 3;
895b8514ebSBarry Smith     if (a[12] == 0.0) {
9077431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",2);
91ede5e988SBarry Smith     }
92ede5e988SBarry Smith 
93ede5e988SBarry Smith     /*
9425783f72SBarry Smith          Now form the inverse
95ede5e988SBarry Smith     */
96ede5e988SBarry Smith 
97ede5e988SBarry Smith    /*     compute inverse(u) */
98ede5e988SBarry Smith 
9925783f72SBarry Smith     for (k = 1; k <= 3; ++k) {
10073d4a2d6SBarry Smith         k3    = 3*k;
10173d4a2d6SBarry Smith         k4    = k3 + k;
10273d4a2d6SBarry Smith 	a[k4] = 1.0 / a[k4];
10339d66777SBarry Smith 	stmp  = -a[k4];
104ede5e988SBarry Smith 	i__2  = k - 1;
10573d4a2d6SBarry Smith         aa    = &a[k3 + 1];
10639d66777SBarry Smith         for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
107ede5e988SBarry Smith 	kp1 = k + 1;
10825783f72SBarry Smith 	if (3 < kp1) continue;
109ede5e988SBarry Smith         ax = aa;
11025783f72SBarry Smith         for (j = kp1; j <= 3; ++j) {
11173d4a2d6SBarry Smith             j3        = 3*j;
11239d66777SBarry Smith 	    stmp      = a[k + j3];
11339d66777SBarry Smith 	    a[k + j3] = 0.0;
11473d4a2d6SBarry Smith             ay        = &a[j3 + 1];
115ede5e988SBarry Smith             for (ll=0; ll<k; ll++) {
11639d66777SBarry Smith               ay[ll] += stmp*ax[ll];
117ede5e988SBarry Smith             }
118ede5e988SBarry Smith 	}
119ede5e988SBarry Smith     }
120ede5e988SBarry Smith 
121ede5e988SBarry Smith    /*    form inverse(u)*inverse(l) */
122ede5e988SBarry Smith 
12325783f72SBarry Smith     for (kb = 1; kb <= 2; ++kb) {
12425783f72SBarry Smith 	k   = 3 - kb;
12573d4a2d6SBarry Smith         k3  = 3*k;
126ede5e988SBarry Smith 	kp1 = k + 1;
12773d4a2d6SBarry Smith         aa  = a + k3;
12825783f72SBarry Smith 	for (i = kp1; i <= 3; ++i) {
129b48ee343SBarry Smith             work[i-1] = aa[i];
13073d4a2d6SBarry Smith 	    aa[i]   = 0.0;
131ede5e988SBarry Smith 	}
13225783f72SBarry Smith 	for (j = kp1; j <= 3; ++j) {
133b48ee343SBarry Smith 	    stmp  = work[j-1];
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 	}
140da10e913SBarry Smith 	l = ipvt[k-1];
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