xref: /petsc/src/mat/impls/baij/seq/dgefa3.c (revision 26fbe8dc1a3c99fb8dddfa572c8c6b3b4ce3ca53)
1be1d678aSKris Buschelman 
2ede5e988SBarry Smith /*
325783f72SBarry Smith      Inverts 3 by 3 matrix using partial pivoting.
471c5468dSBarry Smith 
571c5468dSBarry Smith        Used by the sparse factorization routines in
6dd882469SBarry Smith      src/mat/impls/baij/seq
771c5468dSBarry Smith 
871c5468dSBarry Smith 
971c5468dSBarry Smith        This is a combination of the Linpack routines
1071c5468dSBarry Smith     dgefa() and dgedi() specialized for a size of 3.
1171c5468dSBarry Smith 
12ede5e988SBarry Smith */
13c6db04a5SJed Brown #include <petscsys.h>
14ede5e988SBarry Smith 
154a2ae208SSatish Balay #undef __FUNCT__
1696b95a6bSBarry Smith #define __FUNCT__ "PetscKernel_A_gets_inverse_A_3"
1796b95a6bSBarry Smith PetscErrorCode PetscKernel_A_gets_inverse_A_3(MatScalar *a,PetscReal shift)
18ede5e988SBarry Smith {
19690b6cddSBarry Smith   PetscInt  i__2,i__3,kp1,j,k,l,ll,i,ipvt[3],kb,k3;
20690b6cddSBarry Smith   PetscInt  k4,j3;
21b48ee343SBarry Smith   MatScalar *aa,*ax,*ay,work[9],stmp;
22329f5518SBarry Smith   MatReal   tmp,max;
23ede5e988SBarry Smith 
24ede5e988SBarry Smith /*     gaussian elimination with partial pivoting */
25ede5e988SBarry Smith 
263a40ed3dSBarry Smith   PetscFunctionBegin;
270daa89e9SHong Zhang   shift = .333*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[4]) + PetscAbsScalar(a[8]));
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;
46da10e913SBarry Smith     ipvt[k-1] = l;
47ede5e988SBarry Smith 
48ab055debSShri Abhyankar     if (a[l + k3] == 0.0) {
49*26fbe8dcSKarl Rupp       if (shift == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
50*26fbe8dcSKarl Rupp       else {
51ab055debSShri Abhyankar         /* Shift is applied to single diagonal entry */
52ab055debSShri Abhyankar         a[l + k3] = shift;
53ab055debSShri Abhyankar       }
54ab055debSShri Abhyankar     }
55ede5e988SBarry Smith /*           interchange if necessary */
56ede5e988SBarry Smith 
57ede5e988SBarry Smith     if (l != k) {
5839d66777SBarry Smith       stmp      = a[l + k3];
5973d4a2d6SBarry Smith       a[l + k3] = a[k4];
6039d66777SBarry Smith       a[k4]     = stmp;
61ede5e988SBarry Smith     }
62ede5e988SBarry Smith 
63ede5e988SBarry Smith /*           compute multipliers */
64ede5e988SBarry Smith 
6539d66777SBarry Smith     stmp = -1. / a[k4];
66ede5e988SBarry Smith     i__2 = 3 - k;
6773d4a2d6SBarry Smith     aa   = &a[1 + k4];
68*26fbe8dcSKarl Rupp     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
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];
83*26fbe8dcSKarl Rupp       for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll];
84ede5e988SBarry Smith     }
85ede5e988SBarry Smith   }
86da10e913SBarry Smith   ipvt[2] = 3;
8765e19b50SBarry Smith   if (a[12] == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",2);
88ede5e988SBarry Smith 
89ede5e988SBarry Smith   /*
9025783f72SBarry Smith        Now form the inverse
91ede5e988SBarry Smith   */
92ede5e988SBarry Smith 
93ede5e988SBarry Smith   /*     compute inverse(u) */
94ede5e988SBarry Smith 
9525783f72SBarry Smith   for (k = 1; k <= 3; ++k) {
9673d4a2d6SBarry Smith     k3    = 3*k;
9773d4a2d6SBarry Smith     k4    = k3 + k;
9873d4a2d6SBarry Smith     a[k4] = 1.0 / a[k4];
9939d66777SBarry Smith     stmp  = -a[k4];
100ede5e988SBarry Smith     i__2  = k - 1;
10173d4a2d6SBarry Smith     aa    = &a[k3 + 1];
10239d66777SBarry Smith     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
103ede5e988SBarry Smith     kp1 = k + 1;
10425783f72SBarry Smith     if (3 < kp1) continue;
105ede5e988SBarry Smith     ax = aa;
10625783f72SBarry Smith     for (j = kp1; j <= 3; ++j) {
10773d4a2d6SBarry Smith       j3        = 3*j;
10839d66777SBarry Smith       stmp      = a[k + j3];
10939d66777SBarry Smith       a[k + j3] = 0.0;
11073d4a2d6SBarry Smith       ay        = &a[j3 + 1];
111*26fbe8dcSKarl Rupp       for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll];
112ede5e988SBarry Smith     }
113ede5e988SBarry Smith   }
114ede5e988SBarry Smith 
115ede5e988SBarry Smith   /*    form inverse(u)*inverse(l) */
116ede5e988SBarry Smith 
11725783f72SBarry Smith   for (kb = 1; kb <= 2; ++kb) {
11825783f72SBarry Smith     k   = 3 - kb;
11973d4a2d6SBarry Smith     k3  = 3*k;
120ede5e988SBarry Smith     kp1 = k + 1;
12173d4a2d6SBarry Smith     aa  = a + k3;
12225783f72SBarry Smith     for (i = kp1; i <= 3; ++i) {
123b48ee343SBarry Smith       work[i-1] = aa[i];
12473d4a2d6SBarry Smith       aa[i]     = 0.0;
125ede5e988SBarry Smith     }
12625783f72SBarry Smith     for (j = kp1; j <= 3; ++j) {
127b48ee343SBarry Smith       stmp   = work[j-1];
12825783f72SBarry Smith       ax     = &a[3*j + 1];
12973d4a2d6SBarry Smith       ay     = &a[k3 + 1];
13039d66777SBarry Smith       ay[0] += stmp*ax[0];
13139d66777SBarry Smith       ay[1] += stmp*ax[1];
13239d66777SBarry Smith       ay[2] += stmp*ax[2];
133ede5e988SBarry Smith     }
134da10e913SBarry Smith     l = ipvt[k-1];
135ede5e988SBarry Smith     if (l != k) {
13673d4a2d6SBarry Smith       ax   = &a[k3 + 1];
13725783f72SBarry Smith       ay   = &a[3*l + 1];
13873d4a2d6SBarry Smith       stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
13973d4a2d6SBarry Smith       stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
14073d4a2d6SBarry Smith       stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
141ede5e988SBarry Smith     }
142ede5e988SBarry Smith   }
1433a40ed3dSBarry Smith   PetscFunctionReturn(0);
144ede5e988SBarry Smith }
145ede5e988SBarry Smith 
14671c5468dSBarry Smith 
147