CASM  1.1.0
A Clusters Approach to Statistical Mechanics
ConfigIOHull.cc
Go to the documentation of this file.
2 
3 #include <functional>
4 
6 #include "casm/clex/Norm.hh"
7 #include "casm/clex/PrimClex.hh"
9 
10 namespace CASM {
11 
12 namespace ConfigIO {
13 
17  options["atom_frac"] = Hull::CalculatorPair(
19  options["comp"] =
21  return options;
22 }
23 
27  Clex clex;
28  options["comp"] = Hull::CalculatorPair(Comp().clone(), clex.clone());
29  clex.parse_args("formation_energy,per_species");
30  options["atom_frac"] =
32  return options;
33 }
34 
36 std::string default_hull_calculator() { return "atom_frac"; }
37 
38 // --- BaseHull Template Definitions (only needed here) -------------------
39 
55 template <typename ValueType>
56 BaseHull<ValueType>::BaseHull(const std::string &_name,
57  const std::string &_desc,
58  const std::string &_default_selection,
59  const std::string &_default_composition_type,
60  const Hull::CalculatorOptions &_calculator_map,
61  double _singular_value_tol,
62  double _bottom_facet_tol)
63  : BaseValueFormatter<ValueType, Configuration>(_name, _desc),
64  m_calculator_map(_calculator_map),
65  m_singular_value_tol(_singular_value_tol),
66  m_bottom_facet_tol(_bottom_facet_tol),
67  m_selection(_default_selection),
68  m_composition_type(_default_composition_type),
69  m_initialized(false) {}
70 
78 template <typename ValueType>
79 bool BaseHull<ValueType>::init(const Configuration &_tmplt) const {
80  DB::Selection<Configuration> selection(_tmplt.primclex().db<Configuration>(),
81  m_selection);
82 
83  // Hull typedefs:
84  // typedef std::pair<notstd::cloneable_ptr<CompCalculator>,
85  // notstd::cloneable_ptr<EnergyCalculator> > CalculatorPair;
86  // typedef std::map<std::string, CalculatorPair> CalculatorOptions;
87 
88  auto res = m_calculator_map.find(m_composition_type)->second;
89  res.first->init(_tmplt);
90  res.second->init(_tmplt);
91 
92  m_hull = std::make_shared<Hull>(selection, *res.first, *res.second,
93  m_singular_value_tol, m_bottom_facet_tol);
94  return true;
95 }
96 
104 template <typename ValueType>
106  const Configuration &_config) const {
107  std::stringstream t_ss;
108  t_ss << this->name() << "(" << m_selection << "," << m_composition_type
109  << ")";
110  return t_ss.str();
111 }
112 
135 template <typename ValueType>
136 bool BaseHull<ValueType>::parse_args(const std::string &args) {
137  if (m_initialized) {
138  return false;
139  }
140 
141  std::vector<std::string> splt_vec;
142  boost::split(splt_vec, args, boost::is_any_of(","), boost::token_compress_on);
143 
144  if (splt_vec.size() > 4) {
145  throw std::runtime_error("Attempted to initialize format tag " +
146  this->name() + " with " +
147  std::to_string(splt_vec.size()) + " arguments (" +
148  args + "), but only up to 4 arguments allowed.\n");
149  return false;
150  }
151 
152  while (splt_vec.size() < 4) {
153  splt_vec.push_back("");
154  }
155 
156  m_selection = splt_vec[0].empty() ? m_selection : splt_vec[0];
157  m_composition_type = splt_vec[1].empty() ? m_composition_type : splt_vec[1];
158  m_singular_value_tol =
159  splt_vec[2].empty() ? m_singular_value_tol : std::stod(splt_vec[2]);
160  m_bottom_facet_tol =
161  splt_vec[3].empty() ? m_bottom_facet_tol : std::stod(splt_vec[3]);
162  m_initialized = true;
163 
164  auto it = m_calculator_map.find(m_composition_type);
165  if (it == m_calculator_map.end()) {
166  std::stringstream ss;
167  ss << "Composition type '" << m_composition_type
168  << "' is not recognized. Options are: ";
169  for (auto it = m_calculator_map.begin(); it != m_calculator_map.end();
170  ++it) {
171  std::cerr << it->first << " ";
172  if (it != m_calculator_map.begin()) {
173  ss << ", ";
174  }
175  ss << "'" << it->first << "'";
176  }
177  throw std::runtime_error(ss.str());
178  }
179 
180  // prevent combining formatters
181  return false;
182 }
183 
185 template <typename ValueType>
187  return *m_hull;
188 }
189 
190 // --- OnHull Definitions -------------------
191 
192 const std::string OnHull::Name = "on_hull";
193 
194 const std::string OnHull::Desc =
195  "Whether configuration is a vertex on the formation_energy convex hull "
196  "(i.e., is a groundstate)."
197  " Only one Configuration out of a set that have identical or almost "
198  "identical points in"
199  " composition/energy space will return true."
200  " Accepts arguments ($selection,$composition,$dim_tol,$bottom_tol)."
201  " ($selection may be one of: <filename>, 'ALL', 'CALCULATED', 'MASTER' "
202  "<--default)"
203  " ($composition may be one of: 'comp', 'atom_frac' <--default)"
204  " ($dim_tol: tolerance for detecting composition dimensionality, "
205  "default=1e-8)"
206  " ($bottom_tol: tolerance for detecting which facets form the convex hull "
207  "bottom, default=1e-8)"
208  " For 'comp', 'formation_energy' is used. For 'atom_frac', "
209  "'formation_energy_per_atom' is used."
210  " Ex: on_hull, on_hull(MASTER,comp).";
211 
213 OnHull::OnHull() : BaseHull<bool>(Name, Desc) {}
214 
216 bool OnHull::validate(const Configuration &_config) const {
217  return has_formation_energy(_config);
218 }
219 
225 bool OnHull::evaluate(const Configuration &_config) const {
226  orgQhull::QhullVertexList vertices = _hull().data().vertexList();
227  for (auto it = vertices.begin(); it != vertices.end(); ++it) {
228  if (_hull().configuration(*it).name() == _config.name()) {
229  return true;
230  }
231  }
232  return false;
233 }
234 
235 // --- HullDist Definitions -------------------
236 
237 const std::string HullDist::Name = "hull_dist";
238 
239 const std::string HullDist::Desc =
240  "Distance, in eV, of a configuration's formation_energy_per_atom above the "
241  "convex hull."
242  " Accepts arguments ($selection,$composition,$dim_tol,$bottom_tol)."
243  " ($selection may be one of: <filename>, 'ALL', 'CALCULATED', 'MASTER' "
244  "<--default)"
245  " ($composition may be one of: 'comp', 'atom_frac' <--default)."
246  " ($dim_tol: tolerance for detecting composition dimensionality, "
247  "default=1e-8)"
248  " ($bottom_tol: tolerance for detecting which facets form the convex hull "
249  "bottom, default=1e-8)"
250  " For 'comp', 'formation_energy' is used. For 'atom_frac', "
251  "'formation_energy_per_atom' is used."
252  " Ex: hull_dist, hull_dist(MASTER,comp).";
253 
255 HullDist::HullDist() : BaseHull<double>(Name, Desc) {}
256 
258 bool HullDist::validate(const Configuration &_config) const {
259  return has_formation_energy(_config);
260 }
261 
263 double HullDist::evaluate(const Configuration &_config) const {
264  double d = _hull().dist_to_hull(_config);
265  d = (std::abs(d) < m_dist_to_hull_tol) ? 0.0 : d;
266  return d;
267 }
268 
269 // --- OnClexHull Definitions -------------------
270 
271 const std::string OnClexHull::Name = "on_clex_hull";
272 
273 const std::string OnClexHull::Desc =
274  "Whether configuration is a vertex on the *cluster-expanded* "
275  "formation_energy "
276  "convex hull (i.e., is a *predicted* groundstate)."
277  " Only one Configuration out of a set that have identical or almost "
278  "identical points in"
279  " composition/energy space will return true."
280  " Accepts arguments ($selection,$composition,$dim_tol,$bottom_tol)."
281  " ($selection may be one of: <filename>, 'ALL', 'CALCULATED', 'MASTER' "
282  "<--default)"
283  " ($composition may be one of: 'comp', 'atom_frac' <--default)"
284  " ($dim_tol: tolerance for detecting composition dimensionality, "
285  "default=1e-8)"
286  " ($bottom_tol: tolerance for detecting which facets form the convex hull "
287  "bottom, default=1e-8)"
288  " For 'comp', 'clex(formation_energy)' is used. For 'atom_frac', "
289  "'clex(formation_energy_per_atom)' is used."
290  " Ex: clex_hull_dist, clex_hull_dist(MASTER,comp).";
291 
294  : BaseHull<bool>(Name, Desc, "MASTER", default_hull_calculator(),
296 
301 bool OnClexHull::validate(const Configuration &_config) const { return true; }
302 
308 bool OnClexHull::evaluate(const Configuration &_config) const {
309  orgQhull::QhullVertexList vertices = _hull().data().vertexList();
310  for (auto it = vertices.begin(); it != vertices.end(); ++it) {
311  if (_hull().configuration(*it).name() == _config.name()) {
312  return true;
313  }
314  }
315  return false;
316 }
317 
318 // --- ClexHullDist Definitions -------------------
319 
320 const std::string ClexHullDist::Name = "clex_hull_dist";
321 
322 const std::string ClexHullDist::Desc =
323  "Distance, in eV, of a configuration's *cluster-expanded* "
324  "formation_energy_per_atom above the convex hull."
325  " Accepts arguments ($selection,$composition,$dim_tol,$bottom_tol)."
326  " ($selection may be one of: <filename>, 'ALL', 'CALCULATED', 'MASTER' "
327  "<--default)"
328  " ($composition may be one of: 'comp', 'atom_frac' <--default)"
329  " ($dim_tol: tolerance for detecting composition dimensionality, "
330  "default=1e-8)"
331  " ($bottom_tol: tolerance for detecting which facets form the convex hull "
332  "bottom, default=1e-8)"
333  " For 'comp', 'clex(formation_energy)' is used. For 'atom_frac', "
334  "'clex(formation_energy_per_atom)' is used."
335  " Ex: clex_hull_dist, clex_hull_dist(MASTER,comp).";
336 
339  : BaseHull<double>(Name, Desc, "MASTER", default_hull_calculator(),
341 
346 bool ClexHullDist::validate(const Configuration &_config) const { return true; }
347 
349 double ClexHullDist::evaluate(const Configuration &_config) const {
350  double d = _hull().dist_to_hull(_config);
351  d = (std::abs(d) < m_dist_to_hull_tol) ? 0.0 : d;
352  return d;
353 }
354 
355 template class BaseHull<bool>;
356 template class BaseHull<double>;
357 } // namespace ConfigIO
358 } // namespace CASM
Base class for creating scalar DatumFormatter.
Base class for hull info formatters.
Definition: ConfigIOHull.hh:29
bool parse_args(const std::string &args) override
Determine the selection to use to generate the hull.
constexpr static double m_dist_to_hull_tol
Definition: ConfigIOHull.hh:57
std::string short_header(const Configuration &_config) const override
column header to use
BaseHull(const std::string &_name, const std::string &_desc, const std::string &_default_selection="MASTER", const std::string &_default_composition_type=default_hull_calculator(), const Hull::CalculatorOptions &_calculator_map=hull_calculator_options(), double _singular_value_tol=1e-8, double _bottom_facet_tol=1e-8)
Constructor.
Definition: ConfigIOHull.cc:56
const Hull & _hull() const
const Access the Hull object
bool init(const Configuration &_tmplt) const override
Calculates the convex hull.
Definition: ConfigIOHull.cc:79
Returns predicted formation energy.
Definition: ConfigIO.hh:385
bool parse_args(const std::string &args) override
Expects 'clex', 'clex(formation_energy)', or 'clex(formation_energy_per_species)'.
Definition: ConfigIO.cc:518
std::unique_ptr< Clex > clone() const
Clone using copy constructor.
Definition: ConfigIO.cc:491
virtual double evaluate(const Configuration &_config) const override
Return the distance to the hull.
virtual bool validate(const Configuration &_config) const override
Validate that the Configuration has a cluster expanded formation energy per species.
static const std::string Name
static const std::string Desc
Calculate param composition of a Configuration.
Definition: ConfigIO.hh:72
double evaluate(const Configuration &_config) const override
Return the distance to the hull.
static const std::string Name
static const std::string Desc
bool validate(const Configuration &_config) const override
Validate that the Configuration has a formation energy per species.
virtual bool validate(const Configuration &_config) const override
Validate that the Configuration has a cluster expanded formation energy per species.
virtual bool evaluate(const Configuration &_config) const override
Check if the Configuration is a hull vertex.
static const std::string Desc
static const std::string Name
static const std::string Desc
bool evaluate(const Configuration &_config) const override
Check if the Configuration is a hull vertex.
static const std::string Name
bool validate(const Configuration &_config) const override
Validate that the Configuration has a formation energy per species.
Generate and inspect the convex hull generated from a selection of Configurations.
Definition: Hull.hh:17
std::pair< notstd::cloneable_ptr< CompCalculator >, notstd::cloneable_ptr< EnergyCalculator > > CalculatorPair
Definition: Hull.hh:23
std::map< std::string, CalculatorPair > CalculatorOptions
Definition: Hull.hh:24
double dist_to_hull(const Configuration &config) const
The distance a Configuration is above the hull along the energy axis.
Definition: Hull.cc:231
const orgQhull::Qhull & data() const
const Access the hull object directly
Definition: Hull.cc:119
bool has_formation_energy(const Configuration &_config)
std::string to_string(ENUM val)
Return string representation of enum class.
Definition: io_traits.hh:172
ConfigIO::GenericConfigFormatter< double > formation_energy()
Definition: ConfigIO.cc:670
ConfigIO::GenericConfigFormatter< double > formation_energy_per_species()
Definition: ConfigIO.cc:678
AtomFrac SpeciesFrac
In the future, AtomFrac will actually be atoms only.
Definition: ConfigIO.hh:191
Hull::CalculatorOptions hull_calculator_options()
Returns a map containing default hull calculators.
Definition: ConfigIOHull.cc:15
std::string default_hull_calculator()
Returns "atom_frac".
Definition: ConfigIOHull.cc:36
Hull::CalculatorOptions clex_hull_calculator_options()
Returns a map containing default clex hull calculators.
Definition: ConfigIOHull.cc:25
Main CASM namespace.
Definition: APICommand.hh:8
GenericDatumFormatter< std::string, DataObject > name()
std::unique_ptr< T > clone(const T &obj)