Actual source code: dgefa6.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  2: /*
  3:       Inverts 6 by 6 matrix using gaussian elimination with partial pivoting.

  5:        Used by the sparse factorization routines in
  6:      src/mat/impls/baij/seq

  8:        This is a combination of the Linpack routines
  9:     dgefa() and dgedi() specialized for a size of 6.

 11: */
 12: #include <petscsys.h>

 16: PETSC_EXTERN PetscErrorCode PetscKernel_A_gets_inverse_A_6(MatScalar *a,PetscReal shift,PetscBool allowzeropivot,PetscBool *zeropivotdetected)
 17: {
 18:   PetscInt  i__2,i__3,kp1,j,k,l,ll,i,ipvt[6],kb,k3;
 19:   PetscInt  k4,j3;
 20:   MatScalar *aa,*ax,*ay,work[36],stmp;
 21:   MatReal   tmp,max;

 24:   if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
 25:   shift = .25*shift*(1.e-12 + PetscAbsScalar(a[0]) + PetscAbsScalar(a[7]) + PetscAbsScalar(a[14]) + PetscAbsScalar(a[21]) + PetscAbsScalar(a[28]) + PetscAbsScalar(a[35]));

 27:   /* Parameter adjustments */
 28:   a -= 7;

 30:   for (k = 1; k <= 5; ++k) {
 31:     kp1 = k + 1;
 32:     k3  = 6*k;
 33:     k4  = k3 + k;

 35:     /* find l = pivot index */
 36:     i__2 = 7 - k;
 37:     aa   = &a[k4];
 38:     max  = PetscAbsScalar(aa[0]);
 39:     l    = 1;
 40:     for (ll=1; ll<i__2; ll++) {
 41:       tmp = PetscAbsScalar(aa[ll]);
 42:       if (tmp > max) { max = tmp; l = ll+1;}
 43:     }
 44:     l        += k - 1;
 45:     ipvt[k-1] = l;

 47:     if (a[l + k3] == 0.0) {
 48:       if (shift == 0.0) {
 49:         if (allowzeropivot) {
 51:           PetscInfo1(NULL,"Zero pivot, row %D\n",k-1);
 52:           if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
 53:         } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
 54:       } else {
 55:         /* SHIFT is applied to SINGLE diagonal entry; does this make any sense? */
 56:         a[l + k3] = shift;
 57:       }
 58:     }

 60:     /* interchange if necessary */
 61:     if (l != k) {
 62:       stmp      = a[l + k3];
 63:       a[l + k3] = a[k4];
 64:       a[k4]     = stmp;
 65:     }

 67:     /* compute multipliers */
 68:     stmp = -1. / a[k4];
 69:     i__2 = 6 - k;
 70:     aa   = &a[1 + k4];
 71:     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;

 73:     /* row elimination with column indexing */
 74:     ax = &a[k4+1];
 75:     for (j = kp1; j <= 6; ++j) {
 76:       j3   = 6*j;
 77:       stmp = a[l + j3];
 78:       if (l != k) {
 79:         a[l + j3] = a[k + j3];
 80:         a[k + j3] = stmp;
 81:       }

 83:       i__3 = 6 - k;
 84:       ay   = &a[1+k+j3];
 85:       for (ll=0; ll<i__3; ll++) ay[ll] += stmp*ax[ll];
 86:     }
 87:   }
 88:   ipvt[5] = 6;
 89:   if (a[42] == 0.0) {
 90:     if (allowzeropivot) {
 92:       PetscInfo1(NULL,"Zero pivot, row %D\n",5);
 93:       if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
 94:     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",5);
 95:   }

 97:   /* Now form the inverse */
 98:   /* compute inverse(u) */
 99:   for (k = 1; k <= 6; ++k) {
100:     k3    = 6*k;
101:     k4    = k3 + k;
102:     a[k4] = 1.0 / a[k4];
103:     stmp  = -a[k4];
104:     i__2  = k - 1;
105:     aa    = &a[k3 + 1];
106:     for (ll=0; ll<i__2; ll++) aa[ll] *= stmp;
107:     kp1 = k + 1;
108:     if (6 < kp1) continue;
109:     ax = aa;
110:     for (j = kp1; j <= 6; ++j) {
111:       j3        = 6*j;
112:       stmp      = a[k + j3];
113:       a[k + j3] = 0.0;
114:       ay        = &a[j3 + 1];
115:       for (ll=0; ll<k; ll++) ay[ll] += stmp*ax[ll];
116:     }
117:   }

119:   /* form inverse(u)*inverse(l) */
120:   for (kb = 1; kb <= 5; ++kb) {
121:     k   = 6 - kb;
122:     k3  = 6*k;
123:     kp1 = k + 1;
124:     aa  = a + k3;
125:     for (i = kp1; i <= 6; ++i) {
126:       work[i-1] = aa[i];
127:       aa[i]     = 0.0;
128:     }
129:     for (j = kp1; j <= 6; ++j) {
130:       stmp   = work[j-1];
131:       ax     = &a[6*j + 1];
132:       ay     = &a[k3 + 1];
133:       ay[0] += stmp*ax[0];
134:       ay[1] += stmp*ax[1];
135:       ay[2] += stmp*ax[2];
136:       ay[3] += stmp*ax[3];
137:       ay[4] += stmp*ax[4];
138:       ay[5] += stmp*ax[5];
139:     }
140:     l = ipvt[k-1];
141:     if (l != k) {
142:       ax   = &a[k3 + 1];
143:       ay   = &a[6*l + 1];
144:       stmp = ax[0]; ax[0] = ay[0]; ay[0] = stmp;
145:       stmp = ax[1]; ax[1] = ay[1]; ay[1] = stmp;
146:       stmp = ax[2]; ax[2] = ay[2]; ay[2] = stmp;
147:       stmp = ax[3]; ax[3] = ay[3]; ay[3] = stmp;
148:       stmp = ax[4]; ax[4] = ay[4]; ay[4] = stmp;
149:       stmp = ax[5]; ax[5] = ay[5]; ay[5] = stmp;
150:     }
151:   }
152:   return(0);
153: }