/// Verify that the matrix does not contain elements with abs(element) > 1 and
/// has an acceptable number of non-zero elements.
void SymmetryOperationSymbolParser::verifyMatrix(
    const Kernel::IntMatrix &matrix) {
  for (size_t i = 0; i < matrix.numRows(); ++i) {
    if (!isValidMatrixRow(matrix[i], matrix.numCols())) {
      std::ostringstream error;
      error << "Matrix row " << i << " is invalid. Elements: [" << matrix[i][0]
            << ", " << matrix[i][1] << ", " << matrix[i][2] << "]";
      throw Kernel::Exception::ParseError(error.str(), "", 0);
    }
  }
}
/**
 * Returns the Jones faithful representation of a symmetry operation
 *
 * This method generates a Jones faithful string for the given matrix and
 *vector.
 * The string is generated bases on some rules:
 *
 *  - No spaces:
 *          'x + 1/2' -> 'x+1/2'
 *  - Matrix components occur before vector components:
 *          '1/2+x' -> 'x+1/2'
 *  - No leading '+' signs:
 *          '+x' -> 'x'
 *  - If more than one matrix element is present, they are ordered x, y, z:
 *          'y-x' -> '-x+y'
 *
 * If the matrix is not 3x3, an std::runtime_error exception is thrown.
 *
 * @param matrix
 * @param vector
 * @return
 */
std::string SymmetryOperationSymbolParser::getNormalizedIdentifier(
    const Kernel::IntMatrix &matrix, const V3R &vector) {
  if (matrix.numCols() != 3 || matrix.numRows() != 3) {
    throw std::runtime_error("Matrix is not a 3x3 matrix.");
  }

  std::vector<std::string> symbols{"x", "y", "z"};
  std::vector<std::string> components;

  for (size_t r = 0; r < 3; ++r) {
    std::ostringstream currentComponent;

    for (size_t c = 0; c < 3; ++c) {
      if (matrix[r][c] != 0) {
        if (matrix[r][c] < 0) {
          currentComponent << "-";
        } else {
          if (!currentComponent.str().empty()) {
            currentComponent << "+";
          }
        }

        currentComponent << symbols[c];
      }
    }

    if (vector[r] != 0) {
      if (vector[r] > 0) {
        currentComponent << "+";
      }

      if (vector[r].denominator() != 1) {
        currentComponent << vector[r];
      } else {
        currentComponent << vector[r].numerator();
      }
    }

    components.push_back(currentComponent.str());
  }

  return boost::join(components, ",");
}
/// Returns the order of the symmetry operation based on the matrix. From
/// "Introduction to Crystal Growth and Characterization, Benz and Neumann,
/// Wiley, 2014, p. 51."
size_t
SymmetryOperation::getOrderFromMatrix(const Kernel::IntMatrix &matrix) const {
  int trace = matrix.Trace();
  int determinant = matrix.determinant();

  if (determinant == 1) {
    switch (trace) {
    case 3:
      return 1;
    case 2:
      return 6;
    case 1:
      return 4;
    case 0:
      return 3;
    case -1:
      return 2;
    default:
      break;
    }
  } else if (determinant == -1) {
    switch (trace) {
    case -3:
      return 2;
    case -2:
      return 6;
    case -1:
      return 4;
    case 0:
      return 6;
    case 1:
      return 2;
    default:
      break;
    }
  }

  throw std::runtime_error("There is something wrong with supplied matrix.");
}
/**
 * Returns the symmetry axis for the given matrix
 *
 * According to ITA, 11.2 the axis of a symmetry operation can be determined by
 * solving the Eigenvalue problem \f$Wu = u\f$ for rotations or \f$Wu = -u\f$
 * for rotoinversions. This is implemented using the general real non-symmetric
 * eigen-problem solver provided by the GSL.
 *
 * @param matrix :: Matrix of a SymmetryOperation
 * @return Axis of symmetry element.
 */
V3R SymmetryElementWithAxisGenerator::determineAxis(
    const Kernel::IntMatrix &matrix) const {
  gsl_matrix *eigenMatrix = getGSLMatrix(matrix);
  gsl_matrix *identityMatrix =
      getGSLIdentityMatrix(matrix.numRows(), matrix.numCols());

  gsl_eigen_genv_workspace *eigenWs = gsl_eigen_genv_alloc(matrix.numRows());

  gsl_vector_complex *alpha = gsl_vector_complex_alloc(3);
  gsl_vector *beta = gsl_vector_alloc(3);
  gsl_matrix_complex *eigenVectors = gsl_matrix_complex_alloc(3, 3);

  gsl_eigen_genv(eigenMatrix, identityMatrix, alpha, beta, eigenVectors,
                 eigenWs);
  gsl_eigen_genv_sort(alpha, beta, eigenVectors, GSL_EIGEN_SORT_ABS_DESC);

  double determinant = matrix.determinant();

  Kernel::V3D eigenVector;

  for (size_t i = 0; i < matrix.numCols(); ++i) {
    double eigenValue = GSL_REAL(gsl_complex_div_real(
        gsl_vector_complex_get(alpha, i), gsl_vector_get(beta, i)));

    if (fabs(eigenValue - determinant) < 1e-9) {
      for (size_t j = 0; j < matrix.numRows(); ++j) {
        double element = GSL_REAL(gsl_matrix_complex_get(eigenVectors, j, i));

        eigenVector[j] = element;
      }
    }
  }

  eigenVector *= determinant;

  double sumOfElements = eigenVector.X() + eigenVector.Y() + eigenVector.Z();

  if (sumOfElements < 0) {
    eigenVector *= -1.0;
  }

  gsl_matrix_free(eigenMatrix);
  gsl_matrix_free(identityMatrix);
  gsl_eigen_genv_free(eigenWs);
  gsl_vector_complex_free(alpha);
  gsl_vector_free(beta);
  gsl_matrix_complex_free(eigenVectors);

  double min = 1.0;
  for (size_t i = 0; i < 3; ++i) {
    double absoluteValue = fabs(eigenVector[i]);
    if (absoluteValue != 0.0 &&
        (eigenVector[i] < min && (absoluteValue - fabs(min)) < 1e-9)) {
      min = eigenVector[i];
    }
  }

  V3R axis;
  for (size_t i = 0; i < 3; ++i) {
    axis[i] = static_cast<int>(boost::math::round(eigenVector[i] / min));
  }

  return axis;
}