xref: /petsc/src/mat/impls/baij/seq/dgefa3.c (revision 0daa89e999f16e42c5a19a4a86cef6555d16e163)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
3ede5e988SBarry Smith /*
425783f72SBarry Smith      Inverts 3 by 3 matrix using partial pivoting.
571c5468dSBarry Smith 
671c5468dSBarry Smith        Used by the sparse factorization routines in
7dd882469SBarry Smith      src/mat/impls/baij/seq
871c5468dSBarry Smith 
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 */
14d382aafbSBarry Smith #include "petscsys.h"
15ede5e988SBarry Smith 
164a2ae208SSatish Balay #undef __FUNCT__
174a2ae208SSatish Balay #define __FUNCT__ "Kernel_A_gets_inverse_A_3"
18ab055debSShri Abhyankar PetscErrorCode Kernel_A_gets_inverse_A_3(MatScalar *a,PetscReal shift)
19ede5e988SBarry Smith {
20690b6cddSBarry Smith     PetscInt   i__2,i__3,kp1,j,k,l,ll,i,ipvt[3],kb,k3;
21690b6cddSBarry Smith     PetscInt   k4,j3;
22b48ee343SBarry Smith     MatScalar  *aa,*ax,*ay,work[9],stmp;
23329f5518SBarry Smith     MatReal    tmp,max;
24ede5e988SBarry Smith 
25ede5e988SBarry Smith /*     gaussian elimination with partial pivoting */
26ede5e988SBarry Smith 
273a40ed3dSBarry Smith     PetscFunctionBegin;
28*0daa89e9SHong Zhang     shift = .333*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[4]) + PetscAbsScalar(a[8]));
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 
49ab055debSShri Abhyankar 	if (a[l + k3] == 0.0) {
50ab055debSShri Abhyankar 	  if (shift == 0.0) {
51ab055debSShri Abhyankar 	    SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
52ab055debSShri Abhyankar 	  } else {
53ab055debSShri Abhyankar 	    /* Shift is applied to single diagonal entry */
54ab055debSShri Abhyankar 	    a[l + k3] = shift;
55ab055debSShri Abhyankar 	  }
56ab055debSShri Abhyankar 	}
57ede5e988SBarry Smith /*           interchange if necessary */
58ede5e988SBarry Smith 
59ede5e988SBarry Smith 	if (l != k) {
6039d66777SBarry Smith 	  stmp      = a[l + k3];
6173d4a2d6SBarry Smith 	  a[l + k3] = a[k4];
6239d66777SBarry Smith 	  a[k4]     = stmp;
63ede5e988SBarry Smith         }
64ede5e988SBarry Smith 
65ede5e988SBarry Smith /*           compute multipliers */
66ede5e988SBarry Smith 
6739d66777SBarry Smith 	stmp = -1. / a[k4];
68ede5e988SBarry Smith 	i__2 = 3 - k;
6973d4a2d6SBarry Smith         aa = &a[1 + k4];
70ede5e988SBarry Smith         for (ll=0; ll<i__2; ll++) {
7139d66777SBarry Smith           aa[ll] *= stmp;
72ede5e988SBarry Smith         }
73ede5e988SBarry Smith 
74ede5e988SBarry Smith /*           row elimination with column indexing */
75ede5e988SBarry Smith 
7673d4a2d6SBarry Smith 	ax = &a[k4+1];
77ede5e988SBarry Smith         for (j = kp1; j <= 3; ++j) {
7873d4a2d6SBarry Smith             j3   = 3*j;
7939d66777SBarry Smith 	    stmp = a[l + j3];
80ede5e988SBarry Smith 	    if (l != k) {
8173d4a2d6SBarry Smith 	      a[l + j3] = a[k + j3];
8239d66777SBarry Smith 	      a[k + j3] = stmp;
83ede5e988SBarry Smith             }
84ede5e988SBarry Smith 
85ede5e988SBarry Smith 	    i__3 = 3 - k;
8639d66777SBarry Smith             ay = &a[1+k+j3];
87ede5e988SBarry Smith             for (ll=0; ll<i__3; ll++) {
8839d66777SBarry Smith               ay[ll] += stmp*ax[ll];
89ede5e988SBarry Smith             }
90ede5e988SBarry Smith 	}
91ede5e988SBarry Smith     }
92da10e913SBarry Smith     ipvt[2] = 3;
9365e19b50SBarry Smith     if (a[12] == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",2);
94ede5e988SBarry Smith 
95ede5e988SBarry Smith     /*
9625783f72SBarry Smith          Now form the inverse
97ede5e988SBarry Smith     */
98ede5e988SBarry Smith 
99ede5e988SBarry Smith    /*     compute inverse(u) */
100ede5e988SBarry Smith 
10125783f72SBarry Smith     for (k = 1; k <= 3; ++k) {
10273d4a2d6SBarry Smith         k3    = 3*k;
10373d4a2d6SBarry Smith         k4    = k3 + k;
10473d4a2d6SBarry Smith 	a[k4] = 1.0 / a[k4];
10539d66777SBarry Smith 	stmp  = -a[k4];
106ede5e988SBarry Smith 	i__2  = k - 1;
10773d4a2d6SBarry Smith         aa    = &a[k3 + 1];
10839d66777SBarry Smith         for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
109ede5e988SBarry Smith 	kp1 = k + 1;
11025783f72SBarry Smith 	if (3 < kp1) continue;
111ede5e988SBarry Smith         ax = aa;
11225783f72SBarry Smith         for (j = kp1; j <= 3; ++j) {
11373d4a2d6SBarry Smith             j3        = 3*j;
11439d66777SBarry Smith 	    stmp      = a[k + j3];
11539d66777SBarry Smith 	    a[k + j3] = 0.0;
11673d4a2d6SBarry Smith             ay        = &a[j3 + 1];
117ede5e988SBarry Smith             for (ll=0; ll<k; ll++) {
11839d66777SBarry Smith               ay[ll] += stmp*ax[ll];
119ede5e988SBarry Smith             }
120ede5e988SBarry Smith 	}
121ede5e988SBarry Smith     }
122ede5e988SBarry Smith 
123ede5e988SBarry Smith    /*    form inverse(u)*inverse(l) */
124ede5e988SBarry Smith 
12525783f72SBarry Smith     for (kb = 1; kb <= 2; ++kb) {
12625783f72SBarry Smith 	k   = 3 - kb;
12773d4a2d6SBarry Smith         k3  = 3*k;
128ede5e988SBarry Smith 	kp1 = k + 1;
12973d4a2d6SBarry Smith         aa  = a + k3;
13025783f72SBarry Smith 	for (i = kp1; i <= 3; ++i) {
131b48ee343SBarry Smith             work[i-1] = aa[i];
13273d4a2d6SBarry Smith 	    aa[i]   = 0.0;
133ede5e988SBarry Smith 	}
13425783f72SBarry Smith 	for (j = kp1; j <= 3; ++j) {
135b48ee343SBarry Smith 	    stmp  = work[j-1];
13625783f72SBarry Smith             ax    = &a[3*j + 1];
13773d4a2d6SBarry Smith             ay    = &a[k3 + 1];
13839d66777SBarry Smith             ay[0] += stmp*ax[0];
13939d66777SBarry Smith             ay[1] += stmp*ax[1];
14039d66777SBarry Smith             ay[2] += stmp*ax[2];
141ede5e988SBarry Smith 	}
142da10e913SBarry Smith 	l = ipvt[k-1];
143ede5e988SBarry Smith 	if (l != k) {
14473d4a2d6SBarry Smith             ax = &a[k3 + 1];
14525783f72SBarry Smith             ay = &a[3*l + 1];
14673d4a2d6SBarry Smith             stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
14773d4a2d6SBarry Smith             stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
14873d4a2d6SBarry Smith             stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
149ede5e988SBarry Smith 	}
150ede5e988SBarry Smith     }
1513a40ed3dSBarry Smith     PetscFunctionReturn(0);
152ede5e988SBarry Smith }
153ede5e988SBarry Smith 
15471c5468dSBarry Smith 
155