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