Changeset fed143 in git


Ignore:
Timestamp:
Sep 20, 2010, 10:11:03 AM (14 years ago)
Author:
Frank Seelisch <seelisch@…>
Branches:
(u'spielwiese', 'fe61d9c35bf7c61f2b6cbf1b56e25e2f08d536cc')
Children:
cfd2e6ba5e3f99f46e106336cafaef280361dd61
Parents:
73f78da968d3bf51d9cad24fb4ecf386e3c13276
Message:
rank of a matrix with ground field entries; much faster than linalg.lib's mat_rk

git-svn-id: file:///usr/local/Singular/svn/trunk@13221 2c84dea3-7e68-4137-9b89-c4e89433aadc
Files:
4 edited

Legend:

Unmodified
Added
Removed
  • Singular/iparith.cc

    r73f78d rfed143  
    412412  { "quotient",    0, QUOTIENT_CMD ,      CMD_2},
    413413  { "random",      0, RANDOM_CMD ,        CMD_23},
     414  { "rank",        0, RANK_CMD ,          CMD_12},
    414415  { "read",        0, READ_CMD ,          CMD_12},
    415416  { "reduce",      0, REDUCE_CMD ,        CMD_M},
     
    29602961  res->data =(char *)(long)((i > j) ? i : (rand() % (j-i+1)) + i);
    29612962#endif /* buildin_rand */
     2963  return FALSE;
     2964}
     2965static BOOLEAN jjRANK2(leftv res, leftv u, leftv v)
     2966{
     2967  matrix m =(matrix)u->Data();
     2968  int isRowEchelon = (int)(long)v->Data();
     2969  if (isRowEchelon != 1) isRowEchelon = 0;
     2970  int rank = luRank(m, isRowEchelon);
     2971  res->data =(char *)(long)rank;
    29622972  return FALSE;
    29632973}
     
    39773987,{jjQUOT,      QUOTIENT_CMD,   IDEAL_CMD,      MODUL_CMD,  MODUL_CMD, ALLOW_PLURAL |ALLOW_RING}
    39783988,{jjRANDOM,    RANDOM_CMD,     INT_CMD,        INT_CMD,    INT_CMD, ALLOW_PLURAL |ALLOW_RING}
     3989,{jjRANK2,     RANK_CMD,       INT_CMD,        MATRIX_CMD, INT_CMD, ALLOW_PLURAL |ALLOW_RING}
    39793990,{jjREAD2,     READ_CMD,       STRING_CMD,     LINK_CMD,   STRING_CMD, ALLOW_PLURAL |ALLOW_RING}
    39803991,{jjREDUCE_P,  REDUCE_CMD,     POLY_CMD,       POLY_CMD,   IDEAL_CMD, ALLOW_PLURAL |ALLOW_RING}
     
    50435054  }
    50445055  //res->data = (char *)0;
     5056  return FALSE;
     5057}
     5058static BOOLEAN jjRANK1(leftv res, leftv v)
     5059{
     5060  matrix m =(matrix)v->Data();
     5061  int rank = luRank(m, 0);
     5062  res->data =(char *)(long)rank;
    50455063  return FALSE;
    50465064}
     
    59715989,{kQHWeight,    QHWEIGHT_CMD,    INTVEC_CMD,     MODUL_CMD     , ALLOW_PLURAL |ALLOW_RING}
    59725990,{jjWRONG,      QRING_CMD,       0,              ANY_TYPE      , ALLOW_PLURAL |ALLOW_RING}
     5991,{jjRANK1,      RANK_CMD,        INT_CMD,        MATRIX_CMD    , ALLOW_PLURAL |ALLOW_RING}
    59735992,{jjREAD,       READ_CMD,        STRING_CMD,     LINK_CMD      , ALLOW_PLURAL |ALLOW_RING}
    59745993,{jjREGULARITY, REGULARITY_CMD,  INT_CMD,        LIST_CMD      , NO_PLURAL |ALLOW_RING}
  • Singular/tok.h

    r73f78d rfed143  
    148148  QRING_CMD,
    149149  RANDOM_CMD,
     150  RANK_CMD,
    150151  READ_CMD,
    151152  REPART_CMD,
  • kernel/linearAlgebra.cc

    r73f78d rfed143  
    148148        {
    149149          n = nDiv(pGetCoeff(p), pivotElement);
     150          nNormalize(n);
    150151
    151152          /* filling lMat;
     
    157158          n = nNeg(n);
    158159          for (int cGauss = r + cOffset + 1; cGauss <= cc; cGauss++)
     160          {
    159161            MATELEM(uMat, rGauss, cGauss)
    160162              = pAdd(MATELEM(uMat, rGauss, cGauss),
    161163                     ppMult_nn(MATELEM(uMat, r, cGauss), n));
     164            pNormalize(MATELEM(uMat, rGauss, cGauss));
     165          }
    162166
    163167          nDelete(&n); /* clean up n */
     
    182186bool luInverse(const matrix aMat, matrix &iMat)
    183187{ /* aMat is guaranteed to be an (n x n)-matrix */
    184   int d = aMat->rows();
    185 
    186188  matrix pMat;
    187189  matrix lMat;
     
    196198
    197199  return result;
     200}
     201
     202/* Assumes that aMat is already in row echelon form */
     203int rankFromRowEchelonForm(const matrix aMat)
     204{
     205  int rank = 0;
     206  int rr = aMat->rows(); int cc = aMat->cols();
     207  int r = 1; int c = 1;
     208  while ((r <= rr) && (c <= cc))
     209  {
     210    if (MATELEM(aMat, r, c) == NULL) c++;
     211    else { rank++; r++; }
     212  }
     213  return rank;
     214}
     215
     216int luRank(const matrix aMat, const bool isRowEchelon)
     217{
     218  if (isRowEchelon) return rankFromRowEchelonForm(aMat);
     219  else
     220  { /* compute the LU-decomposition and read off the rank from
     221       the upper triangular matrix of that decomposition */
     222    matrix pMat;
     223    matrix lMat;
     224    matrix uMat;
     225    luDecomp(aMat, pMat, lMat, uMat);
     226    int result = rankFromRowEchelonForm(uMat);
     227 
     228    /* clean-up */
     229    idDelete((ideal*)&pMat);
     230    idDelete((ideal*)&lMat);
     231    idDelete((ideal*)&uMat);
     232 
     233    return result;
     234  }
    198235}
    199236
     
    238275        p = pNeg(p);
    239276        p = pMult(p, pCopy(MATELEM(iMat, r, r)));
     277        pNormalize(p);
    240278        MATELEM(iMat, r, c) = p;
    241279      }
     
    285323        p = pNeg(p);
    286324        p = pMult(p, pCopy(MATELEM(iMat, c, c)));
     325        pNormalize(p);
    287326        MATELEM(iMat, r, c) = p;
    288327      }
     
    349388      p = pAdd(p, ppMult_qq(MATELEM(lMat, r, c), MATELEM(yVec, c, 1) ));
    350389    MATELEM(yVec, r, 1) = pNeg(p);
     390    pNormalize(MATELEM(yVec, r, 1));
    351391  }
    352392 
     
    383423      poly q = pNSet(nInvers(pGetCoeff(MATELEM(uMat, r, nonZeroC))));
    384424      MATELEM(xVec, nonZeroC, 1) = pMult(pNeg(p), q);
     425      pNormalize(MATELEM(xVec, nonZeroC, 1));
    385426      lastNonZeroC = nonZeroC;
    386427    }
  • kernel/linearAlgebra.h

    r73f78d rfed143  
    198198       matrix &iMat       /**< [out] inverse of A if invertible   */
    199199                          );
    200                          
     200
     201/**
     202 * Computes the rank of a given (m x n)-matrix.
     203 *
     204 * The matrix may already be given in row echelon form, e.g., when
     205 * the user has before called luDecomp and passes the upper triangular
     206 * matrix U to luRank. In this case, the second argument can be set to
     207 * true in order to pass this piece of information to the method.
     208 * Otherwise, luDecomp will be called first to compute the matrix U.
     209 * The rank is then read off the matrix U.
     210 *
     211 * @return the rank as an int
     212 * @sa luDecomp(const matrix aMat, matrix &pMat, matrix &lMat, matrix &uMat)
     213 **/
     214int luRank(
     215       const matrix aMat,      /**< [in]  input matrix */
     216       const bool isRowEchelon /**< [in]  if true then aMat is assumed to be
     217                                          already in row echelon form */
     218          );
     219
    201220/**
    202221 * Solves the linear system A*x = b, where A is an (n x n)-matrix
Note: See TracChangeset for help on using the changeset viewer.