xref: /petsc/src/mat/impls/baij/seq/dgefa3.c (revision ce0a2cd1da0658c2b28aad1be2e2c8e41567bece)
1 #define PETSCMAT_DLL
2 
3 /*
4      Inverts 3 by 3 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 3.
13 
14 */
15 #include "petsc.h"
16 
17 #undef __FUNCT__
18 #define __FUNCT__ "Kernel_A_gets_inverse_A_3"
19 PetscErrorCode Kernel_A_gets_inverse_A_3(MatScalar *a,PetscReal Shift)
20 {
21     PetscInt   i__2,i__3,kp1,j,k,l,ll,i,ipvt[3],kb,k3;
22     PetscInt   k4,j3;
23     MatScalar  *aa,*ax,*ay,work[9],stmp;
24     MatReal    tmp,max;
25 
26 /*     gaussian elimination with partial pivoting */
27 
28     PetscFunctionBegin;
29     /* Parameter adjustments */
30     a       -= 4;
31 
32     for (k = 1; k <= 2; ++k) {
33 	kp1 = k + 1;
34         k3  = 3*k;
35         k4  = k3 + k;
36 /*        find l = pivot index */
37 
38 	i__2 = 4 - k;
39         aa = &a[k4];
40         max = PetscAbsScalar(aa[0]);
41         l = 1;
42         for (ll=1; ll<i__2; ll++) {
43           tmp = PetscAbsScalar(aa[ll]);
44           if (tmp > max) { max = tmp; l = ll+1;}
45         }
46         l       += k - 1;
47 	ipvt[k-1] = l;
48 
49 	if (a[l + k3] == 0.0) {
50 	  SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
51 	}
52 
53 /*           interchange if necessary */
54 
55 	if (l != k) {
56 	  stmp      = a[l + k3];
57 	  a[l + k3] = a[k4];
58 	  a[k4]     = stmp;
59         }
60 
61 /*           compute multipliers */
62 
63 	stmp = -1. / a[k4];
64 	i__2 = 3 - k;
65         aa = &a[1 + k4];
66         for (ll=0; ll<i__2; ll++) {
67           aa[ll] *= stmp;
68         }
69 
70 /*           row elimination with column indexing */
71 
72 	ax = &a[k4+1];
73         for (j = kp1; j <= 3; ++j) {
74             j3   = 3*j;
75 	    stmp = a[l + j3];
76 	    if (l != k) {
77 	      a[l + j3] = a[k + j3];
78 	      a[k + j3] = stmp;
79             }
80 
81 	    i__3 = 3 - k;
82             ay = &a[1+k+j3];
83             for (ll=0; ll<i__3; ll++) {
84               ay[ll] += stmp*ax[ll];
85             }
86 	}
87     }
88     ipvt[2] = 3;
89     if (a[12] == 0.0) {
90       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",2);
91     }
92 
93     /*
94          Now form the inverse
95     */
96 
97    /*     compute inverse(u) */
98 
99     for (k = 1; k <= 3; ++k) {
100         k3    = 3*k;
101         k4    = k3 + k;
102 	a[k4] = 1.0 / a[k4];
103 	stmp  = -a[k4];
104 	i__2  = k - 1;
105         aa    = &a[k3 + 1];
106         for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
107 	kp1 = k + 1;
108 	if (3 < kp1) continue;
109         ax = aa;
110         for (j = kp1; j <= 3; ++j) {
111             j3        = 3*j;
112 	    stmp      = a[k + j3];
113 	    a[k + j3] = 0.0;
114             ay        = &a[j3 + 1];
115             for (ll=0; ll<k; ll++) {
116               ay[ll] += stmp*ax[ll];
117             }
118 	}
119     }
120 
121    /*    form inverse(u)*inverse(l) */
122 
123     for (kb = 1; kb <= 2; ++kb) {
124 	k   = 3 - kb;
125         k3  = 3*k;
126 	kp1 = k + 1;
127         aa  = a + k3;
128 	for (i = kp1; i <= 3; ++i) {
129             work[i-1] = aa[i];
130 	    aa[i]   = 0.0;
131 	}
132 	for (j = kp1; j <= 3; ++j) {
133 	    stmp  = work[j-1];
134             ax    = &a[3*j + 1];
135             ay    = &a[k3 + 1];
136             ay[0] += stmp*ax[0];
137             ay[1] += stmp*ax[1];
138             ay[2] += stmp*ax[2];
139 	}
140 	l = ipvt[k-1];
141 	if (l != k) {
142             ax = &a[k3 + 1];
143             ay = &a[3*l + 1];
144             stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
145             stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
146             stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
147 	}
148     }
149     PetscFunctionReturn(0);
150 }
151 
152 
153