xref: /petsc/src/mat/impls/baij/seq/dgefa5.c (revision 4fc747eaadbeca11629f314a99edccbc2ed7b3d3)
1 
2 /*
3       Inverts 5 by 5 matrix using gaussian elimination with 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 5.
10 
11 */
12 #include <petscsys.h>
13 
14 #undef __FUNCT__
15 #define __FUNCT__ "PetscKernel_A_gets_inverse_A_5"
16 PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_5(MatScalar *a,PetscInt *ipvt,MatScalar *work,PetscReal shift,PetscBool allowzeropivot,PetscBool *zeropivotdetected)
17 {
18   PetscInt  i__2,i__3,kp1,j,k,l,ll,i,kb,k3;
19   PetscInt  k4,j3;
20   MatScalar *aa,*ax,*ay,stmp;
21   MatReal   tmp,max;
22 
23   PetscFunctionBegin;
24   if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
25   shift = .25*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[6]) + PetscAbsScalar(a[12]) + PetscAbsScalar(a[18]) + PetscAbsScalar(a[24]));
26 
27   /* Parameter adjustments */
28   a -= 6;
29 
30   for (k = 1; k <= 4; ++k) {
31     kp1 = k + 1;
32     k3  = 5*k;
33     k4  = k3 + k;
34 
35     /* find l = pivot index */
36     i__2 = 6 - 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         if (allowzeropivot) {
50           PetscErrorCode ierr;
51           ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",k-1);CHKERRQ(ierr);
52           *zeropivotdetected = PETSC_TRUE;
53         } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
54       } else {
55         /* SHIFT is applied to SINGLE diagonal entry; does this make any sense? */
56         a[l + k3] = shift;
57       }
58     }
59 
60     /* interchange if necessary */
61     if (l != k) {
62       stmp      = a[l + k3];
63       a[l + k3] = a[k4];
64       a[k4]     = stmp;
65     }
66 
67     /* compute multipliers */
68     stmp = -1. / a[k4];
69     i__2 = 5 - k;
70     aa   = &a[1 + k4];
71     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
72 
73     /* row elimination with column indexing */
74     ax = &a[k4+1];
75     for (j = kp1; j <= 5; ++j) {
76       j3   = 5*j;
77       stmp = a[l + j3];
78       if (l != k) {
79         a[l + j3] = a[k + j3];
80         a[k + j3] = stmp;
81       }
82 
83       i__3 = 5 - k;
84       ay   = &a[1+k+j3];
85       for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll];
86     }
87   }
88   ipvt[4] = 5;
89   if (a[30] == 0.0) {
90     if (allowzeropivot) {
91       PetscErrorCode ierr;
92       ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",4);CHKERRQ(ierr);
93       *zeropivotdetected = PETSC_TRUE;
94     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",4);
95   }
96 
97   /* Now form the inverse */
98   /* compute inverse(u) */
99   for (k = 1; k <= 5; ++k) {
100     k3    = 5*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 (5 < kp1) continue;
109     ax = aa;
110     for (j = kp1; j <= 5; ++j) {
111       j3        = 5*j;
112       stmp      = a[k + j3];
113       a[k + j3] = 0.0;
114       ay        = &a[j3 + 1];
115       for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll];
116     }
117   }
118 
119   /* form inverse(u)*inverse(l) */
120   for (kb = 1; kb <= 4; ++kb) {
121     k   = 5 - kb;
122     k3  = 5*k;
123     kp1 = k + 1;
124     aa  = a + k3;
125     for (i = kp1; i <= 5; ++i) {
126       work[i-1] = aa[i];
127       aa[i]     = 0.0;
128     }
129     for (j = kp1; j <= 5; ++j) {
130       stmp   = work[j-1];
131       ax     = &a[5*j + 1];
132       ay     = &a[k3 + 1];
133       ay[0] += stmp*ax[0];
134       ay[1] += stmp*ax[1];
135       ay[2] += stmp*ax[2];
136       ay[3] += stmp*ax[3];
137       ay[4] += stmp*ax[4];
138     }
139     l = ipvt[k-1];
140     if (l != k) {
141       ax   = &a[k3 + 1];
142       ay   = &a[5*l + 1];
143       stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
144       stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
145       stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
146       stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp;
147       stmp = ax[4]; ax[4] = ay[4]; ay[4] = stmp;
148     }
149   }
150   PetscFunctionReturn(0);
151 }
152 
153