xref: /petsc/src/mat/impls/baij/seq/dgefa5.c (revision 26fbe8dc1a3c99fb8dddfa572c8c6b3b4ce3ca53)
1be1d678aSKris Buschelman 
24224c193SBarry Smith /*
38a36c062SBarry Smith       Inverts 5 by 5 matrix using 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 
144a2ae208SSatish Balay #undef __FUNCT__
1596b95a6bSBarry Smith #define __FUNCT__ "PetscKernel_A_gets_inverse_A_5"
1696b95a6bSBarry Smith PetscErrorCode PetscKernel_A_gets_inverse_A_5(MatScalar *a,PetscInt *ipvt,MatScalar *work,PetscReal shift)
174224c193SBarry Smith {
18766f9fbaSBarry Smith   PetscInt  i__2,i__3,kp1,j,k,l,ll,i,kb,k3;
19690b6cddSBarry Smith   PetscInt  k4,j3;
20766f9fbaSBarry Smith   MatScalar *aa,*ax,*ay,stmp;
21329f5518SBarry Smith   MatReal   tmp,max;
224224c193SBarry Smith 
234224c193SBarry Smith /*     gaussian elimination with partial pivoting */
244224c193SBarry Smith 
253a40ed3dSBarry Smith   PetscFunctionBegin;
26943c8ff5SShri Abhyankar   shift = .25*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[6]) + PetscAbsScalar(a[12]) + PetscAbsScalar(a[18]) + PetscAbsScalar(a[24]));
274224c193SBarry Smith   /* Parameter adjustments */
288a36c062SBarry Smith   a -= 6;
294224c193SBarry Smith 
308a36c062SBarry Smith   for (k = 1; k <= 4; ++k) {
314224c193SBarry Smith     kp1 = k + 1;
328a36c062SBarry Smith     k3  = 5*k;
334224c193SBarry Smith     k4  = k3 + k;
344224c193SBarry Smith /*        find l = pivot index */
354224c193SBarry Smith 
36ed33f8a5SSatish Balay     i__2 = 6 - k;
374224c193SBarry Smith     aa   = &a[k4];
384224c193SBarry Smith     max  = PetscAbsScalar(aa[0]);
394224c193SBarry Smith     l    = 1;
404224c193SBarry Smith     for (ll=1; ll<i__2; ll++) {
414224c193SBarry Smith       tmp = PetscAbsScalar(aa[ll]);
424224c193SBarry Smith       if (tmp > max) { max = tmp; l = ll+1;}
434224c193SBarry Smith     }
444224c193SBarry Smith     l        += k - 1;
45b48ee343SBarry Smith     ipvt[k-1] = l;
464224c193SBarry Smith 
47943c8ff5SShri Abhyankar     if (a[l + k3] == 0.0) {
48*26fbe8dcSKarl Rupp       if (shift == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
49*26fbe8dcSKarl Rupp       else {
50943c8ff5SShri Abhyankar         /* SHIFT is applied to SINGLE diagonal entry; does this make any sense? */
51943c8ff5SShri Abhyankar         a[l + k3] = shift;
52943c8ff5SShri Abhyankar       }
53943c8ff5SShri Abhyankar     }
544224c193SBarry Smith 
554224c193SBarry Smith /*           interchange if necessary */
564224c193SBarry Smith 
574224c193SBarry Smith     if (l != k) {
584224c193SBarry Smith       stmp      = a[l + k3];
594224c193SBarry Smith       a[l + k3] = a[k4];
604224c193SBarry Smith       a[k4]     = stmp;
614224c193SBarry Smith     }
624224c193SBarry Smith 
634224c193SBarry Smith /*           compute multipliers */
644224c193SBarry Smith 
654224c193SBarry Smith     stmp = -1. / a[k4];
668a36c062SBarry Smith     i__2 = 5 - k;
674224c193SBarry Smith     aa   = &a[1 + k4];
68*26fbe8dcSKarl Rupp     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
694224c193SBarry Smith 
704224c193SBarry Smith /*           row elimination with column indexing */
714224c193SBarry Smith 
724224c193SBarry Smith     ax = &a[k4+1];
738a36c062SBarry Smith     for (j = kp1; j <= 5; ++j) {
748a36c062SBarry Smith       j3   = 5*j;
754224c193SBarry Smith       stmp = a[l + j3];
764224c193SBarry Smith       if (l != k) {
774224c193SBarry Smith         a[l + j3] = a[k + j3];
784224c193SBarry Smith         a[k + j3] = stmp;
794224c193SBarry Smith       }
804224c193SBarry Smith 
818a36c062SBarry Smith       i__3 = 5 - k;
824224c193SBarry Smith       ay   = &a[1+k+j3];
83*26fbe8dcSKarl Rupp       for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll];
844224c193SBarry Smith     }
854224c193SBarry Smith   }
86b48ee343SBarry Smith   ipvt[4] = 5;
8765e19b50SBarry Smith   if (a[30] == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",4);
884224c193SBarry Smith 
894224c193SBarry Smith   /*
904224c193SBarry Smith        Now form the inverse
914224c193SBarry Smith   */
924224c193SBarry Smith 
934224c193SBarry Smith   /*     compute inverse(u) */
944224c193SBarry Smith 
958a36c062SBarry Smith   for (k = 1; k <= 5; ++k) {
968a36c062SBarry Smith     k3    = 5*k;
974224c193SBarry Smith     k4    = k3 + k;
984224c193SBarry Smith     a[k4] = 1.0 / a[k4];
994224c193SBarry Smith     stmp  = -a[k4];
1004224c193SBarry Smith     i__2  = k - 1;
1014224c193SBarry Smith     aa    = &a[k3 + 1];
1024224c193SBarry Smith     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
1034224c193SBarry Smith     kp1 = k + 1;
1048a36c062SBarry Smith     if (5 < kp1) continue;
1054224c193SBarry Smith     ax = aa;
1068a36c062SBarry Smith     for (j = kp1; j <= 5; ++j) {
1078a36c062SBarry Smith       j3        = 5*j;
1084224c193SBarry Smith       stmp      = a[k + j3];
1094224c193SBarry Smith       a[k + j3] = 0.0;
1104224c193SBarry Smith       ay        = &a[j3 + 1];
111*26fbe8dcSKarl Rupp       for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll];
1124224c193SBarry Smith     }
1134224c193SBarry Smith   }
1144224c193SBarry Smith 
1154224c193SBarry Smith   /*    form inverse(u)*inverse(l) */
1164224c193SBarry Smith 
1178a36c062SBarry Smith   for (kb = 1; kb <= 4; ++kb) {
1188a36c062SBarry Smith     k   = 5 - kb;
1198a36c062SBarry Smith     k3  = 5*k;
1204224c193SBarry Smith     kp1 = k + 1;
1214224c193SBarry Smith     aa  = a + k3;
1228a36c062SBarry Smith     for (i = kp1; i <= 5; ++i) {
123b48ee343SBarry Smith       work[i-1] = aa[i];
1244224c193SBarry Smith       aa[i]     = 0.0;
1254224c193SBarry Smith     }
1268a36c062SBarry Smith     for (j = kp1; j <= 5; ++j) {
127b48ee343SBarry Smith       stmp   = work[j-1];
1288a36c062SBarry Smith       ax     = &a[5*j + 1];
1294224c193SBarry Smith       ay     = &a[k3 + 1];
1304224c193SBarry Smith       ay[0] += stmp*ax[0];
1314224c193SBarry Smith       ay[1] += stmp*ax[1];
1324224c193SBarry Smith       ay[2] += stmp*ax[2];
1338a36c062SBarry Smith       ay[3] += stmp*ax[3];
1348a36c062SBarry Smith       ay[4] += stmp*ax[4];
1354224c193SBarry Smith     }
136b48ee343SBarry Smith     l = ipvt[k-1];
1374224c193SBarry Smith     if (l != k) {
1384224c193SBarry Smith       ax   = &a[k3 + 1];
1398a36c062SBarry Smith       ay   = &a[5*l + 1];
1404224c193SBarry Smith       stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
1414224c193SBarry Smith       stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
1424224c193SBarry Smith       stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
1438a36c062SBarry Smith       stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp;
1448a36c062SBarry Smith       stmp = ax[4]; ax[4] = ay[4]; ay[4] = stmp;
1454224c193SBarry Smith     }
1464224c193SBarry Smith   }
1473a40ed3dSBarry Smith   PetscFunctionReturn(0);
1484224c193SBarry Smith }
1494224c193SBarry Smith 
150