xref: /petsc/src/mat/impls/baij/seq/dgefa4.c (revision d6acfc2dd824932841bfdf20245c037310c77dc9)
14224c193SBarry Smith /*
2ec1892c8SHong Zhang        Inverts 4 by 4 matrix using gaussian elimination with partial pivoting.
371c5468dSBarry Smith 
471c5468dSBarry Smith        Used by the sparse factorization routines in
5dd882469SBarry Smith      src/mat/impls/baij/seq
671c5468dSBarry Smith 
771c5468dSBarry Smith        This is a combination of the Linpack routines
871c5468dSBarry Smith     dgefa() and dgedi() specialized for a size of 4.
971c5468dSBarry Smith 
104224c193SBarry Smith */
11c6db04a5SJed Brown #include <petscsys.h>
12*d6acfc2dSPierre Jolivet #include <petsc/private/kernels/blockinvert.h>
134224c193SBarry Smith 
14*d6acfc2dSPierre Jolivet PetscErrorCode PetscKernel_A_gets_inverse_A_4(MatScalar *a, PetscReal shift, PetscBool allowzeropivot, PetscBool *zeropivotdetected)
15d71ae5a4SJacob Faibussowitsch {
16690b6cddSBarry Smith   PetscInt   i__2, i__3, kp1, j, k, l, ll, i, ipvt[4], kb, k3;
17690b6cddSBarry Smith   PetscInt   k4, j3;
18b48ee343SBarry Smith   MatScalar *aa, *ax, *ay, work[16], stmp;
19329f5518SBarry Smith   MatReal    tmp, max;
204224c193SBarry Smith 
213a40ed3dSBarry Smith   PetscFunctionBegin;
22c80103daSHong Zhang   if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
230daa89e9SHong Zhang   shift = .25 * shift * (1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[5]) + PetscAbsScalar(a[10]) + PetscAbsScalar(a[15]));
24ec1892c8SHong Zhang 
254224c193SBarry Smith   /* Parameter adjustments */
268a36c062SBarry Smith   a -= 5;
274224c193SBarry Smith 
288a36c062SBarry Smith   for (k = 1; k <= 3; ++k) {
294224c193SBarry Smith     kp1 = k + 1;
308a36c062SBarry Smith     k3  = 4 * k;
314224c193SBarry Smith     k4  = k3 + k;
324224c193SBarry Smith 
33ec1892c8SHong Zhang     /* find l = pivot index */
34ed33f8a5SSatish Balay     i__2 = 5 - k;
354224c193SBarry Smith     aa   = &a[k4];
364224c193SBarry Smith     max  = PetscAbsScalar(aa[0]);
374224c193SBarry Smith     l    = 1;
384224c193SBarry Smith     for (ll = 1; ll < i__2; ll++) {
394224c193SBarry Smith       tmp = PetscAbsScalar(aa[ll]);
409371c9d4SSatish Balay       if (tmp > max) {
419371c9d4SSatish Balay         max = tmp;
429371c9d4SSatish Balay         l   = ll + 1;
439371c9d4SSatish Balay       }
444224c193SBarry Smith     }
454224c193SBarry Smith     l += k - 1;
46da10e913SBarry Smith     ipvt[k - 1] = l;
474224c193SBarry Smith 
485b8514ebSBarry Smith     if (a[l + k3] == 0.0) {
49ec1892c8SHong Zhang       if (shift == 0.0) {
50ec1892c8SHong Zhang         if (allowzeropivot) {
519566063dSJacob Faibussowitsch           PetscCall(PetscInfo(NULL, "Zero pivot, row %" PetscInt_FMT "\n", k - 1));
5270d8d27cSBarry Smith           if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
5398921bdaSJacob Faibussowitsch         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row %" PetscInt_FMT, k - 1);
54ec1892c8SHong Zhang       } else {
551bfe5b3dSBarry Smith         /* SHIFT is applied to SINGLE diagonal entry; does this make any sense? */
561bfe5b3dSBarry Smith         a[l + k3] = shift;
571bfe5b3dSBarry Smith       }
584224c193SBarry Smith     }
594224c193SBarry Smith 
604224c193SBarry Smith     /* interchange if necessary */
614224c193SBarry Smith     if (l != k) {
624224c193SBarry Smith       stmp      = a[l + k3];
634224c193SBarry Smith       a[l + k3] = a[k4];
644224c193SBarry Smith       a[k4]     = stmp;
654224c193SBarry Smith     }
664224c193SBarry Smith 
674224c193SBarry Smith     /* compute multipliers */
684224c193SBarry Smith     stmp = -1. / a[k4];
698a36c062SBarry Smith     i__2 = 4 - k;
704224c193SBarry Smith     aa   = &a[1 + k4];
7126fbe8dcSKarl Rupp     for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp;
724224c193SBarry Smith 
734224c193SBarry Smith     /* row elimination with column indexing */
744224c193SBarry Smith     ax = &a[k4 + 1];
758a36c062SBarry Smith     for (j = kp1; j <= 4; ++j) {
768a36c062SBarry Smith       j3   = 4 * j;
774224c193SBarry Smith       stmp = a[l + j3];
784224c193SBarry Smith       if (l != k) {
794224c193SBarry Smith         a[l + j3] = a[k + j3];
804224c193SBarry Smith         a[k + j3] = stmp;
814224c193SBarry Smith       }
824224c193SBarry Smith 
838a36c062SBarry Smith       i__3 = 4 - k;
844224c193SBarry Smith       ay   = &a[1 + k + j3];
8526fbe8dcSKarl Rupp       for (ll = 0; ll < i__3; ll++) ay[ll] += stmp * ax[ll];
864224c193SBarry Smith     }
874224c193SBarry Smith   }
88da10e913SBarry Smith   ipvt[3] = 4;
892e92ee13SHong Zhang   if (a[20] == 0.0) {
90c0aa6a63SJacob Faibussowitsch     if (PetscLikely(allowzeropivot)) {
919566063dSJacob Faibussowitsch       PetscCall(PetscInfo(NULL, "Zero pivot, row 3\n"));
9270d8d27cSBarry Smith       if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
93c0aa6a63SJacob Faibussowitsch     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row 3");
942e92ee13SHong Zhang   }
954224c193SBarry Smith 
96ec1892c8SHong Zhang   /* Now form the inverse */
974224c193SBarry Smith   /* compute inverse(u) */
988a36c062SBarry Smith   for (k = 1; k <= 4; ++k) {
998a36c062SBarry Smith     k3    = 4 * k;
1004224c193SBarry Smith     k4    = k3 + k;
1014224c193SBarry Smith     a[k4] = 1.0 / a[k4];
1024224c193SBarry Smith     stmp  = -a[k4];
1034224c193SBarry Smith     i__2  = k - 1;
1044224c193SBarry Smith     aa    = &a[k3 + 1];
1054224c193SBarry Smith     for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp;
1064224c193SBarry Smith     kp1 = k + 1;
1078a36c062SBarry Smith     if (4 < kp1) continue;
1084224c193SBarry Smith     ax = aa;
1098a36c062SBarry Smith     for (j = kp1; j <= 4; ++j) {
1108a36c062SBarry Smith       j3        = 4 * j;
1114224c193SBarry Smith       stmp      = a[k + j3];
1124224c193SBarry Smith       a[k + j3] = 0.0;
1134224c193SBarry Smith       ay        = &a[j3 + 1];
11426fbe8dcSKarl Rupp       for (ll = 0; ll < k; ll++) ay[ll] += stmp * ax[ll];
1154224c193SBarry Smith     }
1164224c193SBarry Smith   }
1174224c193SBarry Smith 
1184224c193SBarry Smith   /* form inverse(u)*inverse(l) */
1198a36c062SBarry Smith   for (kb = 1; kb <= 3; ++kb) {
1208a36c062SBarry Smith     k   = 4 - kb;
1218a36c062SBarry Smith     k3  = 4 * k;
1224224c193SBarry Smith     kp1 = k + 1;
1234224c193SBarry Smith     aa  = a + k3;
1248a36c062SBarry Smith     for (i = kp1; i <= 4; ++i) {
125b48ee343SBarry Smith       work[i - 1] = aa[i];
1264224c193SBarry Smith       aa[i]       = 0.0;
1274224c193SBarry Smith     }
1288a36c062SBarry Smith     for (j = kp1; j <= 4; ++j) {
129b48ee343SBarry Smith       stmp = work[j - 1];
1308a36c062SBarry Smith       ax   = &a[4 * j + 1];
1314224c193SBarry Smith       ay   = &a[k3 + 1];
1324224c193SBarry Smith       ay[0] += stmp * ax[0];
1334224c193SBarry Smith       ay[1] += stmp * ax[1];
1344224c193SBarry Smith       ay[2] += stmp * ax[2];
1358a36c062SBarry Smith       ay[3] += stmp * ax[3];
1364224c193SBarry Smith     }
137da10e913SBarry Smith     l = ipvt[k - 1];
1384224c193SBarry Smith     if (l != k) {
1394224c193SBarry Smith       ax    = &a[k3 + 1];
1408a36c062SBarry Smith       ay    = &a[4 * l + 1];
1419371c9d4SSatish Balay       stmp  = ax[0];
1429371c9d4SSatish Balay       ax[0] = ay[0];
1439371c9d4SSatish Balay       ay[0] = stmp;
1449371c9d4SSatish Balay       stmp  = ax[1];
1459371c9d4SSatish Balay       ax[1] = ay[1];
1469371c9d4SSatish Balay       ay[1] = stmp;
1479371c9d4SSatish Balay       stmp  = ax[2];
1489371c9d4SSatish Balay       ax[2] = ay[2];
1499371c9d4SSatish Balay       ay[2] = stmp;
1509371c9d4SSatish Balay       stmp  = ax[3];
1519371c9d4SSatish Balay       ax[3] = ay[3];
1529371c9d4SSatish Balay       ay[3] = stmp;
1534224c193SBarry Smith     }
1544224c193SBarry Smith   }
1553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1564224c193SBarry Smith }
157