7 #include "casm/external/qhull/libqhullcpp/PointCoordinates.h"
8 #include "casm/external/qhull/libqhullcpp/Qhull.h"
21 typedef std::pair<notstd::cloneable_ptr<CompCalculator>,
31 double _singular_value_tol = 1e-14,
double _bottom_facet_tol = 1e-14);
34 const orgQhull::Qhull &
data()
const;
42 const Eigen::Transpose<const Eigen::MatrixXd>
expand()
const;
83 bool operator()(orgQhull::QhullVertex A, orgQhull::QhullVertex B)
const {
84 return A.id() < B.id();
Returns fraction of all species that are a particular species, excluding vacancies.
Generate and inspect the convex hull generated from a selection of Configurations.
const Eigen::Transpose< const Eigen::MatrixXd > expand() const
Orthogonal transformation matrix from a point in the dimension-reduced comp/energy space to the full ...
const Configuration & configuration(const orgQhull::QhullPoint &point) const
Return the configuration corresponding to any point.
std::set< orgQhull::QhullVertex, CompareVertex > m_bottom_vertices
Eigen::VectorXd reduced_point(const Configuration &config) const
Return a vector corresponding to the coordinate of a given configuration in the reduced composition/e...
std::vector< std::pair< orgQhull::QhullFacet, double > > m_bottom_facets
ScalarAttribute< Configuration > EnergyCalculator
const Configuration & groundstate(const Eigen::VectorXd &mu) const
Return the 0K ground state corresponding to the input chemical potential.
DB::Selection< Configuration > m_selection
Hull(const DB::Selection< Configuration > &_selection, const CompCalculator &_comp_calculator=ConfigIO::SpeciesFrac(), const EnergyCalculator &_energy_calculator=ConfigIO::formation_energy_per_species(), double _singular_value_tol=1e-14, double _bottom_facet_tol=1e-14)
Constructor for convex hull in composition/energy space.
std::pair< notstd::cloneable_ptr< CompCalculator >, notstd::cloneable_ptr< EnergyCalculator > > CalculatorPair
const Eigen::MatrixXd & reduce() const
Orthogonal transformation matrix from a point in full comp/energy space to dimension-reduced comp/ene...
notstd::cloneable_ptr< EnergyCalculator > m_energy_calculator
Eigen::VectorXd mu(const orgQhull::QhullFacet facet) const
Return the chemical potential corresponding to a facet.
std::map< std::string, CalculatorPair > CalculatorOptions
Eigen::VectorXd composition(const Configuration &config) const
Use the CompCalculator to return the composition of a Configuration.
VectorXdAttribute< Configuration > CompCalculator
double dist_to_hull(const Configuration &config) const
The distance a Configuration is above the hull along the energy axis.
double energy(const Configuration &config) const
Use the EnergyCalculator to return the energy of a Configuration.
Eigen::VectorXd point(const Configuration &config) const
Return a vector corresponding to the coordinate of a given configuration in full composition/energy s...
notstd::cloneable_ptr< CompCalculator > m_comp_calculator
const orgQhull::Qhull & data() const
const Access the hull object directly
A 'cloneable_ptr' can be used in place of 'unique_ptr'.
ConfigIO::GenericConfigFormatter< jsonParser > config()
ConfigIO::GenericConfigFormatter< double > formation_energy_per_species()
bool operator()(orgQhull::QhullVertex A, orgQhull::QhullVertex B) const