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" 162e92ee13SHong Zhang PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_5(MatScalar *a,PetscInt *ipvt,MatScalar *work,PetscReal shift,PetscBool allowzeropivot,PetscBool *zeropivotdetected) 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; 26*c80103daSHong Zhang if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE; 27*c80103daSHong Zhang 28943c8ff5SShri Abhyankar shift = .25*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[6]) + PetscAbsScalar(a[12]) + PetscAbsScalar(a[18]) + PetscAbsScalar(a[24])); 294224c193SBarry Smith /* Parameter adjustments */ 308a36c062SBarry Smith a -= 6; 314224c193SBarry Smith 328a36c062SBarry Smith for (k = 1; k <= 4; ++k) { 334224c193SBarry Smith kp1 = k + 1; 348a36c062SBarry Smith k3 = 5*k; 354224c193SBarry Smith k4 = k3 + k; 364224c193SBarry Smith /* find l = pivot index */ 374224c193SBarry Smith 38ed33f8a5SSatish Balay i__2 = 6 - k; 394224c193SBarry Smith aa = &a[k4]; 404224c193SBarry Smith max = PetscAbsScalar(aa[0]); 414224c193SBarry Smith l = 1; 424224c193SBarry Smith for (ll=1; ll<i__2; ll++) { 434224c193SBarry Smith tmp = PetscAbsScalar(aa[ll]); 444224c193SBarry Smith if (tmp > max) { max = tmp; l = ll+1;} 454224c193SBarry Smith } 464224c193SBarry Smith l += k - 1; 47b48ee343SBarry Smith ipvt[k-1] = l; 484224c193SBarry Smith 49943c8ff5SShri Abhyankar if (a[l + k3] == 0.0) { 5026fbe8dcSKarl Rupp if (shift == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1); 5126fbe8dcSKarl Rupp else { 52943c8ff5SShri Abhyankar /* SHIFT is applied to SINGLE diagonal entry; does this make any sense? */ 53943c8ff5SShri Abhyankar a[l + k3] = shift; 54943c8ff5SShri Abhyankar } 55943c8ff5SShri Abhyankar } 564224c193SBarry Smith 574224c193SBarry Smith /* interchange if necessary */ 584224c193SBarry Smith 594224c193SBarry Smith if (l != k) { 604224c193SBarry Smith stmp = a[l + k3]; 614224c193SBarry Smith a[l + k3] = a[k4]; 624224c193SBarry Smith a[k4] = stmp; 634224c193SBarry Smith } 644224c193SBarry Smith 654224c193SBarry Smith /* compute multipliers */ 664224c193SBarry Smith 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 744224c193SBarry Smith ax = &a[k4+1]; 758a36c062SBarry Smith for (j = kp1; j <= 5; ++j) { 768a36c062SBarry Smith j3 = 5*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 = 5 - 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 } 88b48ee343SBarry Smith ipvt[4] = 5; 892e92ee13SHong Zhang if (a[30] == 0.0) { 902e92ee13SHong Zhang PetscErrorCode ierr; 912e92ee13SHong Zhang if (allowzeropivot) { 922e92ee13SHong Zhang ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",4);CHKERRQ(ierr); 932e92ee13SHong Zhang *zeropivotdetected = PETSC_TRUE; 942e92ee13SHong Zhang } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",4); 952e92ee13SHong Zhang } 964224c193SBarry Smith 974224c193SBarry Smith /* 984224c193SBarry Smith Now form the inverse 994224c193SBarry Smith */ 1004224c193SBarry Smith 1014224c193SBarry Smith /* compute inverse(u) */ 1024224c193SBarry Smith 1038a36c062SBarry Smith for (k = 1; k <= 5; ++k) { 1048a36c062SBarry Smith k3 = 5*k; 1054224c193SBarry Smith k4 = k3 + k; 1064224c193SBarry Smith a[k4] = 1.0 / a[k4]; 1074224c193SBarry Smith stmp = -a[k4]; 1084224c193SBarry Smith i__2 = k - 1; 1094224c193SBarry Smith aa = &a[k3 + 1]; 1104224c193SBarry Smith for (ll=0; ll<i__2; ll++) aa[ll] *= stmp; 1114224c193SBarry Smith kp1 = k + 1; 1128a36c062SBarry Smith if (5 < kp1) continue; 1134224c193SBarry Smith ax = aa; 1148a36c062SBarry Smith for (j = kp1; j <= 5; ++j) { 1158a36c062SBarry Smith j3 = 5*j; 1164224c193SBarry Smith stmp = a[k + j3]; 1174224c193SBarry Smith a[k + j3] = 0.0; 1184224c193SBarry Smith ay = &a[j3 + 1]; 11926fbe8dcSKarl Rupp for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll]; 1204224c193SBarry Smith } 1214224c193SBarry Smith } 1224224c193SBarry Smith 1234224c193SBarry Smith /* form inverse(u)*inverse(l) */ 1244224c193SBarry Smith 1258a36c062SBarry Smith for (kb = 1; kb <= 4; ++kb) { 1268a36c062SBarry Smith k = 5 - kb; 1278a36c062SBarry Smith k3 = 5*k; 1284224c193SBarry Smith kp1 = k + 1; 1294224c193SBarry Smith aa = a + k3; 1308a36c062SBarry Smith for (i = kp1; i <= 5; ++i) { 131b48ee343SBarry Smith work[i-1] = aa[i]; 1324224c193SBarry Smith aa[i] = 0.0; 1334224c193SBarry Smith } 1348a36c062SBarry Smith for (j = kp1; j <= 5; ++j) { 135b48ee343SBarry Smith stmp = work[j-1]; 1368a36c062SBarry Smith ax = &a[5*j + 1]; 1374224c193SBarry Smith ay = &a[k3 + 1]; 1384224c193SBarry Smith ay[0] += stmp*ax[0]; 1394224c193SBarry Smith ay[1] += stmp*ax[1]; 1404224c193SBarry Smith ay[2] += stmp*ax[2]; 1418a36c062SBarry Smith ay[3] += stmp*ax[3]; 1428a36c062SBarry Smith ay[4] += stmp*ax[4]; 1434224c193SBarry Smith } 144b48ee343SBarry Smith l = ipvt[k-1]; 1454224c193SBarry Smith if (l != k) { 1464224c193SBarry Smith ax = &a[k3 + 1]; 1478a36c062SBarry Smith ay = &a[5*l + 1]; 1484224c193SBarry Smith stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp; 1494224c193SBarry Smith stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp; 1504224c193SBarry Smith stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp; 1518a36c062SBarry Smith stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp; 1528a36c062SBarry Smith stmp = ax[4]; ax[4] = ay[4]; ay[4] = stmp; 1534224c193SBarry Smith } 1544224c193SBarry Smith } 1553a40ed3dSBarry Smith PetscFunctionReturn(0); 1564224c193SBarry Smith } 1574224c193SBarry Smith 158