xref: /petsc/src/mat/impls/baij/seq/dgefa6.c (revision 834855d6effb0d027771461c8e947ee1ce5a1e17)
184643e36SBarry Smith /*
2ec1892c8SHong Zhang       Inverts 6 by 6 matrix using gaussian elimination with partial pivoting.
384643e36SBarry Smith 
484643e36SBarry Smith        Used by the sparse factorization routines in
5dd882469SBarry Smith      src/mat/impls/baij/seq
684643e36SBarry Smith 
784643e36SBarry Smith        This is a combination of the Linpack routines
884643e36SBarry Smith     dgefa() and dgedi() specialized for a size of 6.
984643e36SBarry Smith 
1084643e36SBarry Smith */
11c6db04a5SJed Brown #include <petscsys.h>
12d6acfc2dSPierre Jolivet #include <petsc/private/kernels/blockinvert.h>
1384643e36SBarry Smith 
PetscKernel_A_gets_inverse_A_6(MatScalar * a,PetscReal shift,PetscBool allowzeropivot,PetscBool * zeropivotdetected)14d6acfc2dSPierre Jolivet PetscErrorCode PetscKernel_A_gets_inverse_A_6(MatScalar *a, PetscReal shift, PetscBool allowzeropivot, PetscBool *zeropivotdetected)
15d71ae5a4SJacob Faibussowitsch {
16690b6cddSBarry Smith   PetscInt   i__2, i__3, kp1, j, k, l, ll, i, ipvt[6], kb, k3;
17690b6cddSBarry Smith   PetscInt   k4, j3;
18b48ee343SBarry Smith   MatScalar *aa, *ax, *ay, work[36], stmp;
19329f5518SBarry Smith   MatReal    tmp, max;
2084643e36SBarry Smith 
2184643e36SBarry Smith   PetscFunctionBegin;
22c80103daSHong Zhang   if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
23943c8ff5SShri Abhyankar   shift = .25 * shift * (1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[7]) + PetscAbsScalar(a[14]) + PetscAbsScalar(a[21]) + PetscAbsScalar(a[28]) + PetscAbsScalar(a[35]));
24ec1892c8SHong Zhang 
2584643e36SBarry Smith   /* Parameter adjustments */
2684643e36SBarry Smith   a -= 7;
2784643e36SBarry Smith 
2884643e36SBarry Smith   for (k = 1; k <= 5; ++k) {
2984643e36SBarry Smith     kp1 = k + 1;
3084643e36SBarry Smith     k3  = 6 * k;
3184643e36SBarry Smith     k4  = k3 + k;
3284643e36SBarry Smith 
33ec1892c8SHong Zhang     /* find l = pivot index */
34ed33f8a5SSatish Balay     i__2 = 7 - k;
3584643e36SBarry Smith     aa   = &a[k4];
3684643e36SBarry Smith     max  = PetscAbsScalar(aa[0]);
3784643e36SBarry Smith     l    = 1;
3884643e36SBarry Smith     for (ll = 1; ll < i__2; ll++) {
3984643e36SBarry Smith       tmp = PetscAbsScalar(aa[ll]);
409371c9d4SSatish Balay       if (tmp > max) {
419371c9d4SSatish Balay         max = tmp;
429371c9d4SSatish Balay         l   = ll + 1;
439371c9d4SSatish Balay       }
4484643e36SBarry Smith     }
4584643e36SBarry Smith     l += k - 1;
46b48ee343SBarry Smith     ipvt[k - 1] = l;
4784643e36SBarry Smith 
48943c8ff5SShri Abhyankar     if (a[l + k3] == 0.0) {
49ec1892c8SHong Zhang       if (shift == 0.0) {
50*966bd95aSPierre Jolivet         PetscCheck(allowzeropivot, PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row %" PetscInt_FMT, k - 1);
519566063dSJacob Faibussowitsch         PetscCall(PetscInfo(NULL, "Zero pivot, row %" PetscInt_FMT "\n", k - 1));
5270d8d27cSBarry Smith         if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
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     }
5884643e36SBarry Smith 
5984643e36SBarry Smith     /* interchange if necessary */
6084643e36SBarry Smith     if (l != k) {
6184643e36SBarry Smith       stmp      = a[l + k3];
6284643e36SBarry Smith       a[l + k3] = a[k4];
6384643e36SBarry Smith       a[k4]     = stmp;
6484643e36SBarry Smith     }
6584643e36SBarry Smith 
6684643e36SBarry Smith     /* compute multipliers */
6784643e36SBarry Smith     stmp = -1. / a[k4];
6884643e36SBarry Smith     i__2 = 6 - k;
6984643e36SBarry Smith     aa   = &a[1 + k4];
7026fbe8dcSKarl Rupp     for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp;
7184643e36SBarry Smith 
7284643e36SBarry Smith     /* row elimination with column indexing */
7384643e36SBarry Smith     ax = &a[k4 + 1];
7484643e36SBarry Smith     for (j = kp1; j <= 6; ++j) {
7584643e36SBarry Smith       j3   = 6 * j;
7684643e36SBarry Smith       stmp = a[l + j3];
7784643e36SBarry Smith       if (l != k) {
7884643e36SBarry Smith         a[l + j3] = a[k + j3];
7984643e36SBarry Smith         a[k + j3] = stmp;
8084643e36SBarry Smith       }
8184643e36SBarry Smith 
8284643e36SBarry Smith       i__3 = 6 - k;
8384643e36SBarry Smith       ay   = &a[1 + k + j3];
8426fbe8dcSKarl Rupp       for (ll = 0; ll < i__3; ll++) ay[ll] += stmp * ax[ll];
8584643e36SBarry Smith     }
8684643e36SBarry Smith   }
87b48ee343SBarry Smith   ipvt[5] = 6;
882e92ee13SHong Zhang   if (a[42] == 0.0) {
89*966bd95aSPierre Jolivet     PetscCheck(PetscLikely(allowzeropivot), PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot, row 5");
909566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Zero pivot, row 5\n"));
9170d8d27cSBarry Smith     if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
922e92ee13SHong Zhang   }
9384643e36SBarry Smith 
94ec1892c8SHong Zhang   /* Now form the inverse */
9584643e36SBarry Smith   /* compute inverse(u) */
9684643e36SBarry Smith   for (k = 1; k <= 6; ++k) {
9784643e36SBarry Smith     k3    = 6 * k;
9884643e36SBarry Smith     k4    = k3 + k;
9984643e36SBarry Smith     a[k4] = 1.0 / a[k4];
10084643e36SBarry Smith     stmp  = -a[k4];
10184643e36SBarry Smith     i__2  = k - 1;
10284643e36SBarry Smith     aa    = &a[k3 + 1];
10384643e36SBarry Smith     for (ll = 0; ll < i__2; ll++) aa[ll] *= stmp;
10484643e36SBarry Smith     kp1 = k + 1;
10584643e36SBarry Smith     if (6 < kp1) continue;
10684643e36SBarry Smith     ax = aa;
10784643e36SBarry Smith     for (j = kp1; j <= 6; ++j) {
10884643e36SBarry Smith       j3        = 6 * j;
10984643e36SBarry Smith       stmp      = a[k + j3];
11084643e36SBarry Smith       a[k + j3] = 0.0;
11184643e36SBarry Smith       ay        = &a[j3 + 1];
11226fbe8dcSKarl Rupp       for (ll = 0; ll < k; ll++) ay[ll] += stmp * ax[ll];
11384643e36SBarry Smith     }
11484643e36SBarry Smith   }
11584643e36SBarry Smith 
11684643e36SBarry Smith   /* form inverse(u)*inverse(l) */
11784643e36SBarry Smith   for (kb = 1; kb <= 5; ++kb) {
11884643e36SBarry Smith     k   = 6 - kb;
11984643e36SBarry Smith     k3  = 6 * k;
12084643e36SBarry Smith     kp1 = k + 1;
12184643e36SBarry Smith     aa  = a + k3;
12284643e36SBarry Smith     for (i = kp1; i <= 6; ++i) {
1231f1046a5SBarry Smith       work[i - 1] = aa[i];
12484643e36SBarry Smith       aa[i]       = 0.0;
12584643e36SBarry Smith     }
12684643e36SBarry Smith     for (j = kp1; j <= 6; ++j) {
127b48ee343SBarry Smith       stmp = work[j - 1];
12884643e36SBarry Smith       ax   = &a[6 * j + 1];
12984643e36SBarry Smith       ay   = &a[k3 + 1];
13084643e36SBarry Smith       ay[0] += stmp * ax[0];
13184643e36SBarry Smith       ay[1] += stmp * ax[1];
13284643e36SBarry Smith       ay[2] += stmp * ax[2];
13384643e36SBarry Smith       ay[3] += stmp * ax[3];
13484643e36SBarry Smith       ay[4] += stmp * ax[4];
13584643e36SBarry Smith       ay[5] += stmp * ax[5];
13684643e36SBarry Smith     }
137b48ee343SBarry Smith     l = ipvt[k - 1];
13884643e36SBarry Smith     if (l != k) {
13984643e36SBarry Smith       ax    = &a[k3 + 1];
14084643e36SBarry Smith       ay    = &a[6 * 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;
1539371c9d4SSatish Balay       stmp  = ax[4];
1549371c9d4SSatish Balay       ax[4] = ay[4];
1559371c9d4SSatish Balay       ay[4] = stmp;
1569371c9d4SSatish Balay       stmp  = ax[5];
1579371c9d4SSatish Balay       ax[5] = ay[5];
1589371c9d4SSatish Balay       ay[5] = stmp;
15984643e36SBarry Smith     }
16084643e36SBarry Smith   }
1613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16284643e36SBarry Smith }
163