xref: /petsc/src/mat/impls/baij/seq/dgefa5.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1be1d678aSKris Buschelman 
24224c193SBarry Smith /*
3ec1892c8SHong Zhang       Inverts 5 by 5 matrix using gaussian elimination with partial pivoting.
471c5468dSBarry Smith 
571c5468dSBarry Smith        Used by the sparse factorization routines in
6dd882469SBarry Smith      src/mat/impls/baij/seq
771c5468dSBarry Smith 
871c5468dSBarry Smith        This is a combination of the Linpack routines
971c5468dSBarry Smith     dgefa() and dgedi() specialized for a size of 5.
1071c5468dSBarry Smith 
114224c193SBarry Smith */
12c6db04a5SJed Brown #include <petscsys.h>
134224c193SBarry Smith 
14*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_5(MatScalar *a, PetscInt *ipvt, MatScalar *work, PetscReal shift, PetscBool allowzeropivot, PetscBool *zeropivotdetected) {
15766f9fbaSBarry Smith   PetscInt   i__2, i__3, kp1, j, k, l, ll, i, kb, k3;
16690b6cddSBarry Smith   PetscInt   k4, j3;
17766f9fbaSBarry Smith   MatScalar *aa, *ax, *ay, stmp;
18329f5518SBarry Smith   MatReal    tmp, max;
194224c193SBarry Smith 
203a40ed3dSBarry Smith   PetscFunctionBegin;
21c80103daSHong Zhang   if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
22943c8ff5SShri Abhyankar   shift = .25 * shift * (1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[6]) + PetscAbsScalar(a[12]) + PetscAbsScalar(a[18]) + PetscAbsScalar(a[24]));
23ec1892c8SHong Zhang 
244224c193SBarry Smith   /* Parameter adjustments */
258a36c062SBarry Smith   a -= 6;
264224c193SBarry Smith 
278a36c062SBarry Smith   for (k = 1; k <= 4; ++k) {
284224c193SBarry Smith     kp1 = k + 1;
298a36c062SBarry Smith     k3  = 5 * k;
304224c193SBarry Smith     k4  = k3 + k;
314224c193SBarry Smith 
32ec1892c8SHong Zhang     /* find l = pivot index */
33ed33f8a5SSatish Balay     i__2 = 6 - k;
344224c193SBarry Smith     aa   = &a[k4];
354224c193SBarry Smith     max  = PetscAbsScalar(aa[0]);
364224c193SBarry Smith     l    = 1;
374224c193SBarry Smith     for (ll = 1; ll < i__2; ll++) {
384224c193SBarry Smith       tmp = PetscAbsScalar(aa[ll]);
39*9371c9d4SSatish Balay       if (tmp > max) {
40*9371c9d4SSatish Balay         max = tmp;
41*9371c9d4SSatish Balay         l   = ll + 1;
42*9371c9d4SSatish Balay       }
434224c193SBarry Smith     }
444224c193SBarry Smith     l += k - 1;
45b48ee343SBarry Smith     ipvt[k - 1] = l;
464224c193SBarry Smith 
47943c8ff5SShri Abhyankar     if (a[l + k3] == 0.0) {
48ec1892c8SHong Zhang       if (shift == 0.0) {
49ec1892c8SHong Zhang         if (allowzeropivot) {
509566063dSJacob Faibussowitsch           PetscCall(PetscInfo(NULL, "Zero pivot, row %" PetscInt_FMT "\n", k - 1));
51ec1892c8SHong Zhang           *zeropivotdetected = PETSC_TRUE;
5298921bdaSJacob Faibussowitsch         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row %" PetscInt_FMT, k - 1);
53ec1892c8SHong Zhang       } else {
54943c8ff5SShri Abhyankar         /* SHIFT is applied to SINGLE diagonal entry; does this make any sense? */
55943c8ff5SShri Abhyankar         a[l + k3] = shift;
56943c8ff5SShri Abhyankar       }
57943c8ff5SShri Abhyankar     }
584224c193SBarry Smith 
594224c193SBarry Smith     /* interchange if necessary */
604224c193SBarry Smith     if (l != k) {
614224c193SBarry Smith       stmp      = a[l + k3];
624224c193SBarry Smith       a[l + k3] = a[k4];
634224c193SBarry Smith       a[k4]     = stmp;
644224c193SBarry Smith     }
654224c193SBarry Smith 
664224c193SBarry Smith     /* compute multipliers */
674224c193SBarry Smith     stmp = -1. / a[k4];
688a36c062SBarry Smith     i__2 = 5 - k;
694224c193SBarry Smith     aa   = &a[1 + k4];
7026fbe8dcSKarl Rupp     for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp;
714224c193SBarry Smith 
724224c193SBarry Smith     /* row elimination with column indexing */
734224c193SBarry Smith     ax = &a[k4 + 1];
748a36c062SBarry Smith     for (j = kp1; j <= 5; ++j) {
758a36c062SBarry Smith       j3   = 5 * j;
764224c193SBarry Smith       stmp = a[l + j3];
774224c193SBarry Smith       if (l != k) {
784224c193SBarry Smith         a[l + j3] = a[k + j3];
794224c193SBarry Smith         a[k + j3] = stmp;
804224c193SBarry Smith       }
814224c193SBarry Smith 
828a36c062SBarry Smith       i__3 = 5 - k;
834224c193SBarry Smith       ay   = &a[1 + k + j3];
8426fbe8dcSKarl Rupp       for (ll = 0; ll < i__3; ll++) ay[ll] += stmp * ax[ll];
854224c193SBarry Smith     }
864224c193SBarry Smith   }
87b48ee343SBarry Smith   ipvt[4] = 5;
882e92ee13SHong Zhang   if (a[30] == 0.0) {
89c0aa6a63SJacob Faibussowitsch     if (PetscLikely(allowzeropivot)) {
909566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Zero pivot, row 4\n"));
912e92ee13SHong Zhang       *zeropivotdetected = PETSC_TRUE;
92c0aa6a63SJacob Faibussowitsch     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row 4");
932e92ee13SHong Zhang   }
944224c193SBarry Smith 
95ec1892c8SHong Zhang   /* Now form the inverse */
964224c193SBarry Smith   /* compute inverse(u) */
978a36c062SBarry Smith   for (k = 1; k <= 5; ++k) {
988a36c062SBarry Smith     k3    = 5 * k;
994224c193SBarry Smith     k4    = k3 + k;
1004224c193SBarry Smith     a[k4] = 1.0 / a[k4];
1014224c193SBarry Smith     stmp  = -a[k4];
1024224c193SBarry Smith     i__2  = k - 1;
1034224c193SBarry Smith     aa    = &a[k3 + 1];
1044224c193SBarry Smith     for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp;
1054224c193SBarry Smith     kp1 = k + 1;
1068a36c062SBarry Smith     if (5 < kp1) continue;
1074224c193SBarry Smith     ax = aa;
1088a36c062SBarry Smith     for (j = kp1; j <= 5; ++j) {
1098a36c062SBarry Smith       j3        = 5 * j;
1104224c193SBarry Smith       stmp      = a[k + j3];
1114224c193SBarry Smith       a[k + j3] = 0.0;
1124224c193SBarry Smith       ay        = &a[j3 + 1];
11326fbe8dcSKarl Rupp       for (ll = 0; ll < k; ll++) ay[ll] += stmp * ax[ll];
1144224c193SBarry Smith     }
1154224c193SBarry Smith   }
1164224c193SBarry Smith 
1174224c193SBarry Smith   /* form inverse(u)*inverse(l) */
1188a36c062SBarry Smith   for (kb = 1; kb <= 4; ++kb) {
1198a36c062SBarry Smith     k   = 5 - kb;
1208a36c062SBarry Smith     k3  = 5 * k;
1214224c193SBarry Smith     kp1 = k + 1;
1224224c193SBarry Smith     aa  = a + k3;
1238a36c062SBarry Smith     for (i = kp1; i <= 5; ++i) {
124b48ee343SBarry Smith       work[i - 1] = aa[i];
1254224c193SBarry Smith       aa[i]       = 0.0;
1264224c193SBarry Smith     }
1278a36c062SBarry Smith     for (j = kp1; j <= 5; ++j) {
128b48ee343SBarry Smith       stmp = work[j - 1];
1298a36c062SBarry Smith       ax   = &a[5 * j + 1];
1304224c193SBarry Smith       ay   = &a[k3 + 1];
1314224c193SBarry Smith       ay[0] += stmp * ax[0];
1324224c193SBarry Smith       ay[1] += stmp * ax[1];
1334224c193SBarry Smith       ay[2] += stmp * ax[2];
1348a36c062SBarry Smith       ay[3] += stmp * ax[3];
1358a36c062SBarry Smith       ay[4] += stmp * ax[4];
1364224c193SBarry Smith     }
137b48ee343SBarry Smith     l = ipvt[k - 1];
1384224c193SBarry Smith     if (l != k) {
1394224c193SBarry Smith       ax    = &a[k3 + 1];
1408a36c062SBarry Smith       ay    = &a[5 * l + 1];
141*9371c9d4SSatish Balay       stmp  = ax[0];
142*9371c9d4SSatish Balay       ax[0] = ay[0];
143*9371c9d4SSatish Balay       ay[0] = stmp;
144*9371c9d4SSatish Balay       stmp  = ax[1];
145*9371c9d4SSatish Balay       ax[1] = ay[1];
146*9371c9d4SSatish Balay       ay[1] = stmp;
147*9371c9d4SSatish Balay       stmp  = ax[2];
148*9371c9d4SSatish Balay       ax[2] = ay[2];
149*9371c9d4SSatish Balay       ay[2] = stmp;
150*9371c9d4SSatish Balay       stmp  = ax[3];
151*9371c9d4SSatish Balay       ax[3] = ay[3];
152*9371c9d4SSatish Balay       ay[3] = stmp;
153*9371c9d4SSatish Balay       stmp  = ax[4];
154*9371c9d4SSatish Balay       ax[4] = ay[4];
155*9371c9d4SSatish Balay       ay[4] = stmp;
1564224c193SBarry Smith     }
1574224c193SBarry Smith   }
1583a40ed3dSBarry Smith   PetscFunctionReturn(0);
1594224c193SBarry Smith }
160