/** Assume input matrix contains all distances beween sets A and B.
  * The directed Hausdorff distance from A to B will be the maximum of all
  * minimums in each row.
  * The directed Hausdorff distance from B to A will be the maximum of all
  * minimums in each column.
  * The symmetric Hausdorff distance is the max of the two directed distances.
  * \param m1 The matrix containing distances from A to B.
  * \param hd_ab The directed Hausdorff distance from A to B.
  * \param hd_ba The directed Hausdorff distance from B to A.
  * \return the symmetric Hausdorff distance. 
  */
double Analysis_HausdorffDistance::CalcHausdorffFromMatrix(DataSet_2D const& m1,
                                                           double& hd_ab, double& hd_ba)
{
//  if (m1 == 0) {
//    mprinterr("Internal Error: Analysis_HausdorffDistance::(): Matrix set is null.\n");
//    return -1.0;
//  }
  if (m1.Size() < 1) {
    mprinterr("Error: '%s' is empty.\n", m1.legend());
    return -1.0;
  }
  // Hausdorff distance from A to B. 
  hd_ab = 0.0;
  for (unsigned int row = 0; row != m1.Nrows(); row++)
  {
    double minRow = m1.GetElement(0, row);
    for (unsigned int col = 1; col != m1.Ncols(); col++)
      minRow = std::min( minRow, m1.GetElement(col, row) );
    //mprintf("DEBUG: Min row %6u is %12.4f\n", row, minRow);
    hd_ab = std::max( hd_ab, minRow );
  }
  //mprintf("DEBUG: Hausdorff A to B= %12.4f\n", hd_ab);
  // Hausdorff distance from B to A.
  hd_ba = 0.0;
  for (unsigned int col = 0; col != m1.Ncols(); col++)
  {
    double minCol = m1.GetElement(col, 0);
    for (unsigned int row = 1; row != m1.Nrows(); row++)
      minCol = std::min( minCol, m1.GetElement(col, row) );
    //mprintf("DEBUG: Min col %6u is %12.4f\n", col, minCol);
    hd_ba = std::max( hd_ba, minCol);
  }
  //mprintf("DEBUG: Hausdorff B to A= %12.4f\n", hd_ba);
  // Symmetric Hausdorff distance
  double hd = std::max( hd_ab, hd_ba );
    
  return hd;
}
Exemplo n.º 2
0
/** Get eigenvectors and eigenvalues. They will be stored in descending 
  * order (largest eigenvalue first).
  */
int DataSet_Modes::CalcEigen(DataSet_2D const& mIn, int n_to_calc) {
#ifdef NO_MATHLIB
  mprinterr("Error: modes: Compiled without ARPACK/LAPACK/BLAS routines.\n");
  return 1;
#else
  bool eigenvaluesOnly = false;
  int info = 0;
  if (mIn.MatrixKind() != DataSet_2D::HALF) {
    mprinterr("Error: DataSet_Modes: Eigenvector/value calc only for symmetric matrices.\n");
    return 1;
  }
  // If number to calc is 0, assume we want eigenvalues only
  if (n_to_calc < 1) {
    if (n_to_calc == 0) eigenvaluesOnly = true;
    nmodes_ = (int)mIn.Ncols();
  } else
    nmodes_ = n_to_calc;
  if (nmodes_ > (int)mIn.Ncols()) {
    mprintf("Warning: Specified # of eigenvalues to calc (%i) > matrix dimension (%i).\n",
            nmodes_, mIn.Ncols());
    nmodes_ = mIn.Ncols();
    mprintf("Warning: Only calculating %i eigenvalues.\n", nmodes_);
  }
  if (eigenvaluesOnly)
    mprintf("\tCalculating %i eigenvalues only.\n", nmodes_);
  else
    mprintf("\tCalculating %i eigenvectors and eigenvalues.\n", nmodes_);
  // -----------------------------------------------------------------
  if (nmodes_ == (int)mIn.Ncols()) {
    // Calculate all eigenvalues (and optionally eigenvectors). 
    char jobz = 'V'; // Default: Calc both eigenvectors and eigenvalues
    vecsize_ = mIn.Ncols();
    // Check if only calculating eigenvalues
    if (eigenvaluesOnly) {
      jobz = 'N';
      vecsize_ = 1;
    }
    // Set up space to hold eigenvectors
    if (evectors_ != 0) delete[] evectors_;
    if (!eigenvaluesOnly)
      evectors_ = new double[ nmodes_ * vecsize_ ];
    else
      evectors_ = 0;
    // Set up space to hold eigenvalues
    if (evalues_ != 0) delete[] evalues_;
    evalues_ = new double[ nmodes_ ];
    // Create copy of matrix since call to dspev destroys it
    double* mat = mIn.MatrixArray();
    // Lower triangle; not upper since fortran array layout is inverted w.r.t. C/C++
    char uplo = 'L'; 
    // Allocate temporary workspace
    double* work = new double[ 3 * nmodes_ ];
    // NOTE: The call to dspev is supposed to store eigenvectors in columns. 
    //       However as mentioned above fortran array layout is inverted
    //       w.r.t. C/C++ so eigenvectors end up in rows.
    // NOTE: Eigenvalues/vectors are returned in ascending order.
    dspev_(jobz, uplo, nmodes_, mat, evalues_, evectors_, vecsize_, work, info);
    // If no eigenvectors calcd set vecsize to 0
    if (evectors_==0)
      vecsize_ = 0;
    delete[] work;
    delete[] mat;
    if (info != 0) {
      if (info < 0) {
        mprinterr("Internal Error: from dspev: Argument %i had illegal value.\n", -info);
        mprinterr("Args: %c %c %i matrix %x %x %i work %i\n", jobz, uplo, nmodes_,  
                  evalues_, evectors_, vecsize_, info);
      } else { // info > 0
        mprinterr("Internal Error: from dspev: The algorithm failed to converge.\n");
        mprinterr("%i off-diagonal elements of an intermediate tridiagonal form\n", info);
        mprinterr("did not converge to zero.\n");
      }
      return 1;
    }
  // -----------------------------------------------------------------
  } else {
    // Calculate up to n-1 eigenvalues/vectors using the Implicitly Restarted
    // Arnoldi iteration.
    // FIXME: Eigenvectors obtained with this method appear to have signs
    //        flipped compared to full method - is dot product wrong?
    int nelem = mIn.Ncols(); // Dimension of input matrix (N)
    // Allocate memory to store eigenvectors
    vecsize_ = mIn.Ncols();
    int ncv; // # of columns of the matrix V (evectors_), <= N (mIn.Ncols())
    if (evectors_!=0) delete[] evectors_;
    if (nmodes_*2 <= nelem) 
      ncv = nmodes_*2;
    else 
      ncv = nelem;
    evectors_ = new double[ ncv * nelem ];
    // Temporary storage for eigenvectors to avoid memory overlap in dseupd
    double* eigenvectors = new double[ ncv * nelem ];
    // Allocate memory to store eigenvalues
    if ( evalues_ != 0) delete[] evalues_;
    evalues_ = new double[ nelem ] ; // NOTE: Should this be nmodes?
    // Allocate workspace
    double* workd = new double[ 3 * nelem ];
    int lworkl = ncv * (ncv+8); // NOTE: should this be ncv^2 * (ncv+8)
    double* workl = new double[ lworkl ];
    double* resid = new double[ nelem ];
    // Set parameters for dsaupd (Arnolid)
    int ido = 0; // Reverse comm. flag; 0 = first call
    // The iparam array is used to set parameters for the calc.
    int iparam[11];
    std::fill( iparam, iparam + 11, 0 );
    iparam[0] = 1;   // Method for selecting implicit shifts; 1 = exact
    iparam[2] = 300; // Max # of iterations allowed
    iparam[3] = 1;   // blocksize to be used in the recurrence (code works with only 1).
    iparam[6] = 1;   // Type of eigenproblem being solved; 1: A*x = lambda*x
    double tol = 0;  // Stopping criterion (tolerance); 0 = arpack default 
    char bmat = 'I'; // Type of matrix B that defines semi-inner product; I = identity
    char which[2];   // Which of the Ritz values of OP to compute;
    which[0] = 'L';  // 'LA' = compute the NEV largest eigenvalues
    which[1] = 'A';
    // The ipntr array will hold starting locations in workd and workl arrays
    // for matrices/vectors used by the Lanczos iteration.
    int ipntr[11];
    std::fill( ipntr, ipntr + 11, 0 );
    // Create copy of matrix since it will be modified 
    double* mat = mIn.MatrixArray();
    // LOOP
    bool loop = false;
    do {
      if (loop) {
        // Dot products
        double* target = workd + (ipntr[1] - 1); // -1 since fortran indexing starts at 1
        double* vec    = workd + (ipntr[0] - 1);
        std::fill( target, target + nelem, 0 );
        for(int i = 0; i < nelem; i++) {
          for(int j = i; j < nelem; j++) {
            int ind = nelem * i + j - (i * (i + 1)) / 2;
            target[i] += mat[ind] * vec[j];
            if(i != j)
              target[j] += mat[ind] * vec[i];
          }
        }
      }

      dsaupd_(ido, bmat, nelem, which, nmodes_, tol, resid,
              ncv, eigenvectors, nelem, iparam, ipntr, workd, workl,
              lworkl, info);
      loop = (ido == -1 || ido == 1);
    } while ( loop ); // END LOOP

    if (info != 0) {
      mprinterr("Error: DataSet_Modes: dsaupd returned %i\n",info);
    } else {
      int rvec = 1;
      char howmny = 'A';
      double sigma = 0.0;
      int* select = new int[ ncv ];
      dseupd_(rvec, howmny, select, evalues_, evectors_, nelem, sigma,
              bmat, nelem, which, nmodes_, tol, resid,
              ncv, eigenvectors, nelem, iparam, ipntr, workd, workl,
              lworkl, info);
      delete[] select;
    } 
    delete[] mat;
    delete[] workl;
    delete[] workd;
    delete[] resid;
    delete[] eigenvectors;
    if (info != 0) { 
      mprinterr("Error: DataSet_Modes: dseupd returned %i\n",info);
      return 1;
    }
  }
  // Eigenvalues and eigenvectors are in ascending order. Resort so that
  // they are in descending order (i.e. largest eigenvalue first).
  int pivot = nmodes_ / 2;
  int nmode = nmodes_ - 1;
  double* vtmp = 0;
  if (evectors_ != 0) 
    vtmp = new double[ vecsize_ ];
  for (int mode = 0; mode < pivot; ++mode) {
    // Swap eigenvalue
    double eval = evalues_[mode];
    evalues_[mode] = evalues_[nmode];
    evalues_[nmode] = eval;
    // Swap eigenvector
    if (vtmp != 0) {
      double* Vec0 = evectors_ + (mode  * vecsize_);
      double* Vec1 = evectors_ + (nmode * vecsize_);
      std::copy( Vec0, Vec0 + vecsize_, vtmp );
      std::copy( Vec1, Vec1 + vecsize_, Vec0 );
      std::copy( vtmp, vtmp + vecsize_, Vec1 );
    }
    --nmode;
  }
  if (vtmp != 0) delete[] vtmp;

  return 0;
#endif
}