/*! solve overdetermined or underdetermined A*X=Y using dgels
  with the sum of residual squares output\n
  The residual is set as the columnwise sum of residual squares 
  for overdetermined problems
  while it is always zero for underdetermined problems.
*/
inline long dgematrix::dgels(dgematrix& mat, drovector& residual)
{
#ifdef  CPPL_VERBOSE
  std::cerr << "# [MARK] dgematrix::dgels(dgematrix&, drovector&)"
            << std::endl;
#endif//CPPL_VERBOSE
  
#ifdef  CPPL_DEBUG
  if(M!=mat.M){
    std::cerr << "[ERROR] dgematrix::dgels(dgematrix&, drovector&) "
              << std::endl
              << "These two matrices cannot be solved." << std::endl
              << "Your input was (" << M << "x" << N << ") and ("
              << mat.M << "x" << mat.N << ")." << std::endl;
    exit(1);
  }
#endif//CPPL_DEBUG
  
  residual.resize(mat.N); residual.zero();
  
  if(M<N){ //underdetermined
    dgematrix tmp(N,mat.N);
    for(long i=0; i<mat.M; i++){ for(long j=0; j<mat.N; j++){
      tmp(i,j) =mat(i,j);
    }}
    mat.clear();
    swap(mat,tmp);
  }
  
  char TRANS('N');
  long NRHS(mat.N), LDA(M), LDB(mat.M),
    LWORK(min(M,N)+max(min(M,N),NRHS)), INFO(1);
  double *WORK(new double[LWORK]);
  dgels_(TRANS, M, N, NRHS, Array, LDA, mat.Array, LDB, WORK, LWORK, INFO);
  delete [] WORK;
  
  if(M>N){ //overdetermined
    for(long i=0; i<residual.L; i++){ for(long j=0; j<M-N; j++){
      residual(i) += std::pow(mat(N+j,i), 2.0);
    }}
    
    dgematrix tmp(N,mat.N);
    for(long i=0; i<tmp.M; i++){ for(long j=0; j<tmp.N; j++){
      tmp(i,j) =mat(i,j);
    }}
    mat.clear();
    swap(mat,tmp);
  }
  
  if(INFO!=0){
    std::cerr << "[WARNING] dgematrix::dgels(dgematrix&, drovector&) "
              << "Serious trouble happend. INFO = " << INFO << "."
              << std::endl;
  }
  return INFO;
}
/*! calculate the least-squares-least-norm solution for overdetermined or 
  underdetermined A*x=y using dgelss\n
*/
inline long dgematrix::dgelss(dgematrix& B, dcovector& S, long& RANK,
                              const double RCOND =-1. )
{
#ifdef  CPPL_VERBOSE
  std::cerr << "# [MARK] dgematrix::dgelss(dgematrix&, dcovector&, long& const double)"
            << std::endl;
#endif//CPPL_VERBOSE
  
#ifdef  CPPL_DEBUG
  if(M!=B.M){
    std::cerr << "[ERROR] dgematrix::dgelss"
              << "(dgematrix&, dcovector&, long&, const double) " << std::endl
              << "These matrix and vector cannot be solved." << std::endl
              << "Your input was (" << M << "x" << N << ") and ("
              << B.M << "x" << B.N  << ")." << std::endl;
    exit(1);
  }
#endif//CPPL_DEBUG    
  
  if(M<N){ //underdetermined
    dgematrix tmp(N,B.N);
    for(long i=0; i<B.M; i++){
      for(long j=0; j<B.N; j++){
        tmp(i,j)=B(i,j);
      }
    }
    B.clear();
    swap(B,tmp);
  }
  
  S.resize(min(M,N));
  
  long NRHS(B.N), LDA(M), LDB(B.M),
    LWORK(3*min(M,N)+max(max(2*min(M,N),max(M,N)), NRHS)), INFO(1);
  double *WORK(new double[LWORK]);
  dgelss_(M, N, NRHS, Array, LDA, B.Array, LDB, S.Array,
          RCOND, RANK, WORK, LWORK, INFO);
  delete [] WORK;

  if(INFO!=0){
    std::cerr << "[WARNING] dgematrix::dgelss"
              << "(dgematrix&, docvector&, long, const double) "
              << "Serious trouble happend. INFO = " << INFO << "."
              << std::endl;
  }
  return INFO;
}