IPP Software Navigation Tools IPP Links Communication Pan-STARRS Links

Changeset 23983


Ignore:
Timestamp:
Apr 28, 2009, 10:59:22 AM (17 years ago)
Author:
eugene
Message:

replace LU-decomp solver with real Gauss-Jordan elimination (which actually inverts A) -- based on Kahan description : see Ohana/src/libohana/doc/kahan-gji.pdf

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/psLib/src/math/psMatrix.c

    r19843 r23983  
    318318// This used to be "a temporary gauss-jordan solver based on gene's version based on the Numerical Recipes
    319319// version".  However, it's been removed due to copyright, and replaced with LU Decomposition solving.
     320# if (0)
    320321bool psMatrixGJSolve(psImage *a,
    321322                     psVector *b
     
    328329    PS_ASSERT_INT_EQUAL(a->numCols, a->numRows, false);
    329330    PS_ASSERT_VECTOR_SIZE(b, (long int)a->numCols, false);
    330 
    331 # if (0)
    332     // XXX expand this out explicitly below
    333     // Check for non-finite entries
    334 #define MATRIX_CHECK_NONFINITE_CASE(TYPE,MATRIX) \
    335   case PS_TYPE_##TYPE: { \
    336       ps##TYPE **values = MATRIX->data.TYPE; /* Dereference */ \
    337       int numCols = (MATRIX)->numCols, numRows = (MATRIX)->numRows; /* Size of matrix */ \
    338       for (int i = 0; i < numRows; i++) { \
    339           for (int j = 0; j < numCols; j++) { \
    340               if (!isfinite(values[i][j])) { \
    341                   // psError(PS_ERR_BAD_PARAMETER_VALUE, 3,             
    342                   // "Input matrix contains non-finite elements: matrix[%d][%d] is %.2f\n",
    343                   // i, j, values[i][j]);                               
    344                   return false; \
    345               } \
    346           } \
    347       } \
    348       break; \
    349   }
    350 # endif
    351331
    352332    // Check for non-finite entries in matrix
     
    410390}
    411391
     392# else
     393
     394// Gauss-Jordan elimination using full pivots based on Press et al's description.  Substantially
     395// reworked for psLib style: major modifications to conform to C indexing, use a boolean to track the
     396// completed pivot rows and catch the singular matrix early on.  Also, much cleaner control loops
     397// than their implementation.  XXX this really needs to check on round-off errors (see version by
     398// William Kahan
     399bool psMatrixGJSolve(psImage *a,
     400                     psVector *b
     401    )
     402{
     403    PS_ASSERT_IMAGE_NON_NULL(a, false);
     404    PS_ASSERT_VECTOR_NON_NULL(b, false);
     405    PS_ASSERT_IMAGE_TYPE_F32_OR_F64(a, false);
     406    PS_ASSERT_VECTOR_TYPE_F32_OR_F64(b, false);
     407    psAssert (a->type.type == b->type.type, "types must match for now");
     408    PS_ASSERT_INT_EQUAL(a->numCols, a->numRows, false);
     409    PS_ASSERT_VECTOR_SIZE(b, (long int)a->numCols, false);
     410
     411    // Check for non-finite entries in matrix
     412    switch (a->type.type) {
     413      case PS_TYPE_F32: {
     414          psF32 **values = a->data.F32; /* Dereference */
     415          int numCols = a->numCols, numRows = a->numRows; /* Size of matrix */
     416          for (int i = 0; i < numRows; i++) {
     417              for (int j = 0; j < numCols; j++) {
     418                  if (!isfinite(values[i][j])) {
     419                      return false;
     420                  }
     421              }
     422          }
     423          break;
     424      }
     425      case PS_TYPE_F64: {
     426          psF64 **values = a->data.F64; /* Dereference */
     427          int numCols = a->numCols, numRows = a->numRows; /* Size of matrix */
     428          for (int i = 0; i < numRows; i++) {
     429              for (int j = 0; j < numCols; j++) {
     430                  if (!isfinite(values[i][j])) {
     431                      return false;
     432                  }
     433              }
     434          }
     435          break;
     436      }
     437      default:
     438        psAbort("Should never get here.");
     439    }
     440 
     441    // Following the algorithm laid out by Press et al., we loop along the matrix diagonal,
     442    // but we do not operate on the diagonal elements in order instead, we are looking for
     443    // the current max element and operating on that diagonal element.  this is effectively
     444    // column pivoting.  row pivoting is perfomed explicitly. 
     445
     446    int nSquare = a->numCols;
     447
     448    psVector *colIndexV = psVectorAlloc(nSquare, PS_TYPE_S32);
     449    psVector *rowIndexV = psVectorAlloc(nSquare, PS_TYPE_S32);
     450    psVector *pivotV    = psVectorAlloc(nSquare, PS_TYPE_S32);
     451    psVectorInit (pivotV, 0.0);
     452
     453    if (a->type.type == PS_TYPE_F32) {
     454        psF32 **A = a->data.F32;
     455        psF32  *B = b->data.F32;
     456        int *colIndex = colIndexV->data.S32;
     457        int *rowIndex = rowIndexV->data.S32;
     458        int *pivot    = pivotV->data.S32;
     459
     460        for (int diag = 0; diag < nSquare; diag++) {
     461
     462            psF32 maxval = 0.0;
     463            int maxrow = 0;
     464            int maxcol = 0;
     465
     466            // search for the next pivot
     467            for (int row = 0; row < nSquare; row++) {
     468                if (!finite(A[row][diag])) goto escape;
     469
     470                // if we have already operated on this row (pivot[row] is true), skip it
     471                if (pivot[row]) continue;
     472
     473                // if we have not yet operated on this row (pivot[row] is false), look for pivot for this row
     474                for (int col = 0; col < nSquare; col++) {
     475                    if (pivot[col]) continue;
     476                    if (fabs (A[row][col]) < maxval) continue;
     477                    maxval = fabs (A[row][col]);
     478                    maxrow = row;
     479                    maxcol = col;
     480                }
     481            }
     482
     483            // if pivot[maxcol] is set, we have already done this row: this implies a singular matrix
     484            if (pivot[maxcol]) goto escape;
     485            pivot[maxcol] = 1;
     486
     487            // if the selected pivot is off the diagonal, do a row swap
     488            if (maxrow != maxcol) {
     489                for (int col = 0; col < nSquare; col++) PS_SWAP (A[maxrow][col], A[maxcol][col]);
     490                PS_SWAP (B[maxrow], B[maxcol]);
     491            }
     492            rowIndex[diag] = maxrow;
     493            colIndex[diag] = maxcol;
     494            if (A[maxcol][maxcol] == 0.0) goto escape;
     495            // XXX Kahan replaces the 0.0 pivot with epsilon*(largest element in column) + underFlow
     496
     497            /* rescale by pivot reciprocal */
     498            psF32 tmpval = 1.0 / A[maxcol][maxcol];
     499            A[maxcol][maxcol] = 1.0;
     500            for (int col = 0; col < nSquare; col++) A[maxcol][col] *= tmpval;
     501            B[maxcol] *= tmpval;
     502            // XXX measure the pivot growth and trigger on over/under flow
     503            // growth *= tmpval;
     504            // fprintf (stderr, "column: %d, growth: %e, epsilon: %e\n", maxcol, growth, epsilon);
     505
     506            /* adjust the elements above the pivot */
     507            for (int row = 0; row < nSquare; row++) {
     508                if (row == maxcol) continue;
     509                tmpval = A[row][maxcol];
     510                A[row][maxcol] = 0.0;
     511                for (int col = 0; col < nSquare; col++) A[row][col] -= A[maxcol][col]*tmpval;
     512                B[row] -= B[maxcol]*tmpval;
     513            }
     514        }
     515
     516        // swap back the inverse matrix based on the row swaps above
     517        for (int col = nSquare - 1; col >= 0; col--) {
     518            if (rowIndex[col] != colIndex[col]) {
     519                for (int row = 0; row < nSquare; row++) PS_SWAP (A[row][rowIndex[col]], A[row][colIndex[col]]);
     520            }
     521        }
     522    } else {
     523        psF64 **A = a->data.F64;
     524        psF64  *B = b->data.F64;
     525        int *colIndex = colIndexV->data.S32;
     526        int *rowIndex = rowIndexV->data.S32;
     527        int *pivot    = pivotV->data.S32;
     528
     529        for (int diag = 0; diag < nSquare; diag++) {
     530
     531            psF64 maxval = 0.0;
     532            int maxrow = 0;
     533            int maxcol = 0;
     534
     535            // search for the next pivot
     536            for (int row = 0; row < nSquare; row++) {
     537                if (!finite(A[row][diag])) goto escape;
     538
     539                // if we have already operated on this row (pivot[row] is true), skip it
     540                if (pivot[row]) continue;
     541
     542                // if we have not yet operated on this row (pivot[row] is false), look for pivot for this row
     543                for (int col = 0; col < nSquare; col++) {
     544                    if (pivot[col]) continue;
     545                    if (fabs (A[row][col]) < maxval) continue;
     546                    maxval = fabs (A[row][col]);
     547                    maxrow = row;
     548                    maxcol = col;
     549                }
     550            }
     551
     552            // if pivot[maxcol] is set, we have already done this row: this implies a singular matrix
     553            if (pivot[maxcol]) goto escape;
     554            pivot[maxcol] = 1;
     555
     556            // if the selected pivot is off the diagonal, do a row swap
     557            if (maxrow != maxcol) {
     558                for (int col = 0; col < nSquare; col++) PS_SWAP (A[maxrow][col], A[maxcol][col]);
     559                PS_SWAP (B[maxrow], B[maxcol]);
     560            }
     561            rowIndex[diag] = maxrow;
     562            colIndex[diag] = maxcol;
     563            if (A[maxcol][maxcol] == 0.0) goto escape;
     564            // XXX Kahan replaces the 0.0 pivot with epsilon*(largest element in column) + underFlow
     565
     566            /* rescale by pivot reciprocal */
     567            psF64 tmpval = 1.0 / A[maxcol][maxcol];
     568            A[maxcol][maxcol] = 1.0;
     569            for (int col = 0; col < nSquare; col++) A[maxcol][col] *= tmpval;
     570            B[maxcol] *= tmpval;
     571            // XXX measure the pivot growth and trigger on over/under flow
     572            // growth *= tmpval;
     573            // fprintf (stderr, "column: %d, growth: %e, epsilon: %e\n", maxcol, growth, epsilon);
     574
     575            /* adjust the elements above the pivot */
     576            for (int row = 0; row < nSquare; row++) {
     577                if (row == maxcol) continue;
     578                tmpval = A[row][maxcol];
     579                A[row][maxcol] = 0.0;
     580                for (int col = 0; col < nSquare; col++) A[row][col] -= A[maxcol][col]*tmpval;
     581                B[row] -= B[maxcol]*tmpval;
     582            }
     583        }
     584
     585        // swap back the inverse matrix based on the row swaps above
     586        for (int col = nSquare - 1; col >= 0; col--) {
     587            if (rowIndex[col] != colIndex[col]) {
     588                for (int row = 0; row < nSquare; row++) PS_SWAP (A[row][rowIndex[col]], A[row][colIndex[col]]);
     589            }
     590        }
     591    }
     592
     593    psFree (pivotV);
     594    psFree (rowIndexV);
     595    psFree (colIndexV);
     596    return true;
     597
     598escape:
     599    psFree (pivotV);
     600    psFree (rowIndexV);
     601    psFree (colIndexV);
     602    return false;
     603}
     604
     605# endif
    412606
    413607psImage* psMatrixInvert(psImage* out,
Note: See TracChangeset for help on using the changeset viewer.