xref: /petsc/src/mat/impls/baij/seq/dgefa5.c (revision 2e92ee13a8395f820cc1e3fd74a7607ed52efa2a)
1 
2 /*
3       Inverts 5 by 5 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 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   /* gaussian elimination with partial pivoting */
24 
25   PetscFunctionBegin;
26   shift = .25*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[6]) + PetscAbsScalar(a[12]) + PetscAbsScalar(a[18]) + PetscAbsScalar(a[24]));
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     /* find l = pivot index */
35 
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) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
49       else {
50         /* SHIFT is applied to SINGLE diagonal entry; does this make any sense? */
51         a[l + k3] = shift;
52       }
53     }
54 
55     /* interchange if necessary */
56 
57     if (l != k) {
58       stmp      = a[l + k3];
59       a[l + k3] = a[k4];
60       a[k4]     = stmp;
61     }
62 
63     /* compute multipliers */
64 
65     stmp = -1. / a[k4];
66     i__2 = 5 - k;
67     aa   = &a[1 + k4];
68     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
69 
70     /* row elimination with column indexing */
71 
72     ax = &a[k4+1];
73     for (j = kp1; j <= 5; ++j) {
74       j3   = 5*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 = 5 - k;
82       ay   = &a[1+k+j3];
83       for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll];
84     }
85   }
86   ipvt[4] = 5;
87   if (a[30] == 0.0) {
88     PetscErrorCode ierr;
89     if (allowzeropivot) {
90       ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",4);CHKERRQ(ierr);
91       *zeropivotdetected = PETSC_TRUE;
92     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",4);
93   }
94 
95   /*
96    Now form the inverse
97   */
98 
99   /* compute inverse(u) */
100 
101   for (k = 1; k <= 5; ++k) {
102     k3    = 5*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 (5 < kp1) continue;
111     ax = aa;
112     for (j = kp1; j <= 5; ++j) {
113       j3        = 5*j;
114       stmp      = a[k + j3];
115       a[k + j3] = 0.0;
116       ay        = &a[j3 + 1];
117       for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll];
118     }
119   }
120 
121   /* form inverse(u)*inverse(l) */
122 
123   for (kb = 1; kb <= 4; ++kb) {
124     k   = 5 - kb;
125     k3  = 5*k;
126     kp1 = k + 1;
127     aa  = a + k3;
128     for (i = kp1; i <= 5; ++i) {
129       work[i-1] = aa[i];
130       aa[i]     = 0.0;
131     }
132     for (j = kp1; j <= 5; ++j) {
133       stmp   = work[j-1];
134       ax     = &a[5*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       ay[3] += stmp*ax[3];
140       ay[4] += stmp*ax[4];
141     }
142     l = ipvt[k-1];
143     if (l != k) {
144       ax   = &a[k3 + 1];
145       ay   = &a[5*l + 1];
146       stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
147       stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
148       stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
149       stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp;
150       stmp = ax[4]; ax[4] = ay[4]; ay[4] = stmp;
151     }
152   }
153   PetscFunctionReturn(0);
154 }
155 
156