xref: /petsc/src/mat/impls/baij/seq/dgefa6.c (revision 98921bda46e76d7aaed9e0138c5ff9d0ce93f355)
1be1d678aSKris Buschelman 
284643e36SBarry Smith /*
3ec1892c8SHong Zhang       Inverts 6 by 6 matrix using gaussian elimination with partial pivoting.
484643e36SBarry Smith 
584643e36SBarry Smith        Used by the sparse factorization routines in
6dd882469SBarry Smith      src/mat/impls/baij/seq
784643e36SBarry Smith 
884643e36SBarry Smith        This is a combination of the Linpack routines
984643e36SBarry Smith     dgefa() and dgedi() specialized for a size of 6.
1084643e36SBarry Smith 
1184643e36SBarry Smith */
12c6db04a5SJed Brown #include <petscsys.h>
1384643e36SBarry Smith 
142e92ee13SHong Zhang PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_6(MatScalar *a,PetscReal shift,PetscBool allowzeropivot,PetscBool *zeropivotdetected)
1584643e36SBarry Smith {
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]);
4084643e36SBarry Smith       if (tmp > max) { max = tmp; l = ll+1;}
4184643e36SBarry Smith     }
4284643e36SBarry Smith     l        += k - 1;
43b48ee343SBarry Smith     ipvt[k-1] = l;
4484643e36SBarry Smith 
45943c8ff5SShri Abhyankar     if (a[l + k3] == 0.0) {
46ec1892c8SHong Zhang       if (shift == 0.0) {
47ec1892c8SHong Zhang         if (allowzeropivot) {
48ec1892c8SHong Zhang           PetscErrorCode ierr;
49c0aa6a63SJacob Faibussowitsch           ierr = PetscInfo1(NULL,"Zero pivot, row %" PetscInt_FMT "\n",k-1);CHKERRQ(ierr);
5070d8d27cSBarry Smith           if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
51*98921bdaSJacob Faibussowitsch         } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %" PetscInt_FMT,k-1);
52ec1892c8SHong Zhang       } else {
53943c8ff5SShri Abhyankar         /* SHIFT is applied to SINGLE diagonal entry; does this make any sense? */
54943c8ff5SShri Abhyankar         a[l + k3] = shift;
55943c8ff5SShri Abhyankar       }
56943c8ff5SShri Abhyankar     }
5784643e36SBarry Smith 
5884643e36SBarry Smith     /* interchange if necessary */
5984643e36SBarry Smith     if (l != k) {
6084643e36SBarry Smith       stmp      = a[l + k3];
6184643e36SBarry Smith       a[l + k3] = a[k4];
6284643e36SBarry Smith       a[k4]     = stmp;
6384643e36SBarry Smith     }
6484643e36SBarry Smith 
6584643e36SBarry Smith     /* compute multipliers */
6684643e36SBarry Smith     stmp = -1. / a[k4];
6784643e36SBarry Smith     i__2 = 6 - k;
6884643e36SBarry Smith     aa   = &a[1 + k4];
6926fbe8dcSKarl Rupp     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
7084643e36SBarry Smith 
7184643e36SBarry Smith     /* row elimination with column indexing */
7284643e36SBarry Smith     ax = &a[k4+1];
7384643e36SBarry Smith     for (j = kp1; j <= 6; ++j) {
7484643e36SBarry Smith       j3   = 6*j;
7584643e36SBarry Smith       stmp = a[l + j3];
7684643e36SBarry Smith       if (l != k) {
7784643e36SBarry Smith         a[l + j3] = a[k + j3];
7884643e36SBarry Smith         a[k + j3] = stmp;
7984643e36SBarry Smith       }
8084643e36SBarry Smith 
8184643e36SBarry Smith       i__3 = 6 - k;
8284643e36SBarry Smith       ay   = &a[1+k+j3];
8326fbe8dcSKarl Rupp       for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll];
8484643e36SBarry Smith     }
8584643e36SBarry Smith   }
86b48ee343SBarry Smith   ipvt[5] = 6;
872e92ee13SHong Zhang   if (a[42] == 0.0) {
88c0aa6a63SJacob Faibussowitsch     if (PetscLikely(allowzeropivot)) {
89ec1892c8SHong Zhang       PetscErrorCode ierr;
90c0aa6a63SJacob Faibussowitsch       ierr = PetscInfo(NULL,"Zero pivot, row 5\n");CHKERRQ(ierr);
9170d8d27cSBarry Smith       if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
92c0aa6a63SJacob Faibussowitsch     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row 5");
932e92ee13SHong Zhang   }
9484643e36SBarry Smith 
95ec1892c8SHong Zhang   /* Now form the inverse */
9684643e36SBarry Smith   /* compute inverse(u) */
9784643e36SBarry Smith   for (k = 1; k <= 6; ++k) {
9884643e36SBarry Smith     k3    = 6*k;
9984643e36SBarry Smith     k4    = k3 + k;
10084643e36SBarry Smith     a[k4] = 1.0 / a[k4];
10184643e36SBarry Smith     stmp  = -a[k4];
10284643e36SBarry Smith     i__2  = k - 1;
10384643e36SBarry Smith     aa    = &a[k3 + 1];
10484643e36SBarry Smith     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
10584643e36SBarry Smith     kp1 = k + 1;
10684643e36SBarry Smith     if (6 < kp1) continue;
10784643e36SBarry Smith     ax = aa;
10884643e36SBarry Smith     for (j = kp1; j <= 6; ++j) {
10984643e36SBarry Smith       j3        = 6*j;
11084643e36SBarry Smith       stmp      = a[k + j3];
11184643e36SBarry Smith       a[k + j3] = 0.0;
11284643e36SBarry Smith       ay        = &a[j3 + 1];
11326fbe8dcSKarl Rupp       for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll];
11484643e36SBarry Smith     }
11584643e36SBarry Smith   }
11684643e36SBarry Smith 
11784643e36SBarry Smith   /* form inverse(u)*inverse(l) */
11884643e36SBarry Smith   for (kb = 1; kb <= 5; ++kb) {
11984643e36SBarry Smith     k   = 6 - kb;
12084643e36SBarry Smith     k3  = 6*k;
12184643e36SBarry Smith     kp1 = k + 1;
12284643e36SBarry Smith     aa  = a + k3;
12384643e36SBarry Smith     for (i = kp1; i <= 6; ++i) {
1241f1046a5SBarry Smith       work[i-1] = aa[i];
12584643e36SBarry Smith       aa[i]     = 0.0;
12684643e36SBarry Smith     }
12784643e36SBarry Smith     for (j = kp1; j <= 6; ++j) {
128b48ee343SBarry Smith       stmp   = work[j-1];
12984643e36SBarry Smith       ax     = &a[6*j + 1];
13084643e36SBarry Smith       ay     = &a[k3 + 1];
13184643e36SBarry Smith       ay[0] += stmp*ax[0];
13284643e36SBarry Smith       ay[1] += stmp*ax[1];
13384643e36SBarry Smith       ay[2] += stmp*ax[2];
13484643e36SBarry Smith       ay[3] += stmp*ax[3];
13584643e36SBarry Smith       ay[4] += stmp*ax[4];
13684643e36SBarry Smith       ay[5] += stmp*ax[5];
13784643e36SBarry Smith     }
138b48ee343SBarry Smith     l = ipvt[k-1];
13984643e36SBarry Smith     if (l != k) {
14084643e36SBarry Smith       ax   = &a[k3 + 1];
14184643e36SBarry Smith       ay   = &a[6*l + 1];
14284643e36SBarry Smith       stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
14384643e36SBarry Smith       stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
14484643e36SBarry Smith       stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
14584643e36SBarry Smith       stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp;
14684643e36SBarry Smith       stmp = ax[4]; ax[4] = ay[4]; ay[4] = stmp;
14784643e36SBarry Smith       stmp = ax[5]; ax[5] = ay[5]; ay[5] = stmp;
14884643e36SBarry Smith     }
14984643e36SBarry Smith   }
15084643e36SBarry Smith   PetscFunctionReturn(0);
15184643e36SBarry Smith }
15284643e36SBarry Smith 
153