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