1be1d678aSKris Buschelman 284643e36SBarry Smith /* 384643e36SBarry Smith Inverts 7 by 7 matrix using 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 7. 1084643e36SBarry Smith 1184643e36SBarry Smith */ 12c6db04a5SJed Brown #include <petscsys.h> 1384643e36SBarry Smith 144a2ae208SSatish Balay #undef __FUNCT__ 1596b95a6bSBarry Smith #define __FUNCT__ "PetscKernel_A_gets_inverse_A_7" 162e92ee13SHong Zhang PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_7(MatScalar *a,PetscReal shift,PetscBool allowzeropivot,PetscBool *zeropivotdetected) 1784643e36SBarry Smith { 18690b6cddSBarry Smith PetscInt i__2,i__3,kp1,j,k,l,ll,i,ipvt[7],kb,k3; 19690b6cddSBarry Smith PetscInt k4,j3; 201f1046a5SBarry Smith MatScalar *aa,*ax,*ay,work[49],stmp; 21329f5518SBarry Smith MatReal tmp,max; 2284643e36SBarry Smith 2384643e36SBarry Smith /* gaussian elimination with partial pivoting */ 2484643e36SBarry Smith 2584643e36SBarry 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[8]) + PetscAbsScalar(a[16]) + PetscAbsScalar(a[24]) + PetscAbsScalar(a[32]) + PetscAbsScalar(a[40]) + PetscAbsScalar(a[48])); 29943c8ff5SShri Abhyankar 3084643e36SBarry Smith /* Parameter adjustments */ 3184643e36SBarry Smith a -= 8; 3284643e36SBarry Smith 3384643e36SBarry Smith for (k = 1; k <= 6; ++k) { 3484643e36SBarry Smith kp1 = k + 1; 3584643e36SBarry Smith k3 = 7*k; 3684643e36SBarry Smith k4 = k3 + k; 3784643e36SBarry Smith /* find l = pivot index */ 3884643e36SBarry Smith 39ed33f8a5SSatish Balay i__2 = 8 - k; 4084643e36SBarry Smith aa = &a[k4]; 4184643e36SBarry Smith max = PetscAbsScalar(aa[0]); 4284643e36SBarry Smith l = 1; 4384643e36SBarry Smith for (ll=1; ll<i__2; ll++) { 4484643e36SBarry Smith tmp = PetscAbsScalar(aa[ll]); 4584643e36SBarry Smith if (tmp > max) { max = tmp; l = ll+1;} 4684643e36SBarry Smith } 4784643e36SBarry Smith l += k - 1; 481f1046a5SBarry Smith ipvt[k-1] = l; 4984643e36SBarry Smith 50943c8ff5SShri Abhyankar if (a[l + k3] == 0.0) { 5126fbe8dcSKarl Rupp if (shift == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1); 5226fbe8dcSKarl Rupp 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 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 6884643e36SBarry Smith stmp = -1. / a[k4]; 6984643e36SBarry Smith i__2 = 7 - k; 7084643e36SBarry Smith aa = &a[1 + k4]; 7126fbe8dcSKarl Rupp for (ll=0; ll<i__2; ll++) aa[ll] *= stmp; 7284643e36SBarry Smith 7384643e36SBarry Smith /* row elimination with column indexing */ 7484643e36SBarry Smith 7584643e36SBarry Smith ax = &a[k4+1]; 7684643e36SBarry Smith for (j = kp1; j <= 7; ++j) { 7784643e36SBarry Smith j3 = 7*j; 7884643e36SBarry Smith stmp = a[l + j3]; 7984643e36SBarry Smith if (l != k) { 8084643e36SBarry Smith a[l + j3] = a[k + j3]; 8184643e36SBarry Smith a[k + j3] = stmp; 8284643e36SBarry Smith } 8384643e36SBarry Smith 8484643e36SBarry Smith i__3 = 7 - k; 8584643e36SBarry Smith ay = &a[1+k+j3]; 8626fbe8dcSKarl Rupp for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll]; 8784643e36SBarry Smith } 8884643e36SBarry Smith } 891f1046a5SBarry Smith ipvt[6] = 7; 902e92ee13SHong Zhang if (a[56] == 0.0) { 912e92ee13SHong Zhang PetscErrorCode ierr; 922e92ee13SHong Zhang if (allowzeropivot) { 932e92ee13SHong Zhang ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",6);CHKERRQ(ierr); 942e92ee13SHong Zhang *zeropivotdetected = PETSC_TRUE; 952e92ee13SHong Zhang } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",6); 962e92ee13SHong Zhang } 9784643e36SBarry Smith 9884643e36SBarry Smith /* 9984643e36SBarry Smith Now form the inverse 10084643e36SBarry Smith */ 10184643e36SBarry Smith 10284643e36SBarry Smith /* compute inverse(u) */ 10384643e36SBarry Smith 10484643e36SBarry Smith for (k = 1; k <= 7; ++k) { 10584643e36SBarry Smith k3 = 7*k; 10684643e36SBarry Smith k4 = k3 + k; 10784643e36SBarry Smith a[k4] = 1.0 / a[k4]; 10884643e36SBarry Smith stmp = -a[k4]; 10984643e36SBarry Smith i__2 = k - 1; 11084643e36SBarry Smith aa = &a[k3 + 1]; 11184643e36SBarry Smith for (ll=0; ll<i__2; ll++) aa[ll] *= stmp; 11284643e36SBarry Smith kp1 = k + 1; 11384643e36SBarry Smith if (7 < kp1) continue; 11484643e36SBarry Smith ax = aa; 11584643e36SBarry Smith for (j = kp1; j <= 7; ++j) { 11684643e36SBarry Smith j3 = 7*j; 11784643e36SBarry Smith stmp = a[k + j3]; 11884643e36SBarry Smith a[k + j3] = 0.0; 11984643e36SBarry Smith ay = &a[j3 + 1]; 12026fbe8dcSKarl Rupp for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll]; 12184643e36SBarry Smith } 12284643e36SBarry Smith } 12384643e36SBarry Smith 12484643e36SBarry Smith /* form inverse(u)*inverse(l) */ 12584643e36SBarry Smith 12684643e36SBarry Smith for (kb = 1; kb <= 6; ++kb) { 12784643e36SBarry Smith k = 7 - kb; 12884643e36SBarry Smith k3 = 7*k; 12984643e36SBarry Smith kp1 = k + 1; 13084643e36SBarry Smith aa = a + k3; 13184643e36SBarry Smith for (i = kp1; i <= 7; ++i) { 1321f1046a5SBarry Smith work[i-1] = aa[i]; 13384643e36SBarry Smith aa[i] = 0.0; 13484643e36SBarry Smith } 13584643e36SBarry Smith for (j = kp1; j <= 7; ++j) { 1361f1046a5SBarry Smith stmp = work[j-1]; 13784643e36SBarry Smith ax = &a[7*j + 1]; 13884643e36SBarry Smith ay = &a[k3 + 1]; 13984643e36SBarry Smith ay[0] += stmp*ax[0]; 14084643e36SBarry Smith ay[1] += stmp*ax[1]; 14184643e36SBarry Smith ay[2] += stmp*ax[2]; 14284643e36SBarry Smith ay[3] += stmp*ax[3]; 14384643e36SBarry Smith ay[4] += stmp*ax[4]; 14484643e36SBarry Smith ay[5] += stmp*ax[5]; 14584643e36SBarry Smith ay[6] += stmp*ax[6]; 14684643e36SBarry Smith } 1471f1046a5SBarry Smith l = ipvt[k-1]; 14884643e36SBarry Smith if (l != k) { 14984643e36SBarry Smith ax = &a[k3 + 1]; 15015091d37SBarry Smith ay = &a[7*l + 1]; 15184643e36SBarry Smith stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp; 15284643e36SBarry Smith stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp; 15384643e36SBarry Smith stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp; 15484643e36SBarry Smith stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp; 15584643e36SBarry Smith stmp = ax[4]; ax[4] = ay[4]; ay[4] = stmp; 15684643e36SBarry Smith stmp = ax[5]; ax[5] = ay[5]; ay[5] = stmp; 15784643e36SBarry Smith stmp = ax[6]; ax[6] = ay[6]; ay[6] = stmp; 15884643e36SBarry Smith } 15984643e36SBarry Smith } 16084643e36SBarry Smith PetscFunctionReturn(0); 16184643e36SBarry Smith } 16284643e36SBarry Smith 16384643e36SBarry Smith 164