CASM  1.1.0
A Clusters Approach to Statistical Mechanics
CompositionConverter.hh
Go to the documentation of this file.
1 #ifndef CASM_CompositionConverter_HH
2 #define CASM_CompositionConverter_HH
3 
4 #include <map>
5 #include <vector>
6 
8 #include "casm/global/eigen.hh"
9 
10 namespace CASM {
17  public:
18  typedef unsigned int size_type;
19 
22 
24  template <typename ComponentIterator, typename... EndMembers>
25  CompositionConverter(ComponentIterator begin, ComponentIterator end,
26  Eigen::VectorXd _origin, EndMembers... _end_members);
27 
29  template <typename ComponentIterator>
30  CompositionConverter(ComponentIterator begin, ComponentIterator end,
31  Eigen::VectorXd _origin, Eigen::MatrixXd _end_members);
32 
35 
37  static std::string comp_var(size_type i);
38 
40  std::vector<std::string> components() const;
41 
43  Eigen::VectorXd origin() const;
44 
47 
50 
53 
56 
60 
63 
67 
70 
73 
75  std::string mol_formula() const;
76 
78  std::string param_formula() const;
79 
81  std::string origin_formula() const;
82 
84  std::string end_member_formula(size_type i) const;
85 
87  std::string comp_formula(size_type i) const;
88 
91  std::string comp_n_formula(size_type i) const;
92 
95  std::string param_chem_pot_formula(size_type i) const;
96 
97  private:
99  void _add_end_member(Eigen::VectorXd _end_member);
100 
102  template <typename... EndMembers>
103  void _add_end_member(Eigen::VectorXd _end_member, EndMembers... _others);
104 
107  void _check_size(const Eigen::MatrixXd &vec) const;
108 
111 
113  std::string _n_formula(const Eigen::VectorXd &vec) const;
114 
118  std::vector<std::string> m_components;
119 
124 
130 
134 
138 };
139 
142 template <typename OutputIterator>
143 OutputIterator standard_composition_axes(
144  ParamComposition::AllowedOccupants _allowed_occs, OutputIterator result);
145 
148  std::ostream &stream,
149  const std::map<std::string, CompositionConverter> &map);
150 
152 void display_comp(std::ostream &stream, const CompositionConverter &f,
153  int indent = 0);
154 
156 void display_comp_n(std::ostream &stream, const CompositionConverter &f,
157  int indent = 0);
158 
160 void display_param_chem_pot(std::ostream &stream, const CompositionConverter &f,
161  int indent = 0);
162 
167  const ParamComposition::AllowedOccupants &_allowed_occs);
168 
171  const ParamComposition::AllowedOccupants &_allowed_occs,
172  double tol = 1e-14);
173 
177  const ParamComposition::AllowedOccupants &_allowed_occs,
178  double tol = 1e-14);
179 
180 // ------ Definitions ---------------------------------------------
181 
193 template <typename ComponentIterator, typename... EndMembers>
195  ComponentIterator end,
196  Eigen::VectorXd _origin,
197  EndMembers... _end_members)
198  : m_components(begin, end), m_origin(_origin) {
199  _add_end_member(_end_members...);
200 
201  _check_size(_origin);
202 
204 }
205 
218 template <typename ComponentIterator>
220  ComponentIterator end,
221  Eigen::VectorXd _origin,
222  Eigen::MatrixXd _end_members)
223  : m_components(begin, end), m_origin(_origin), m_end_members(_end_members) {
224  _check_size(_origin);
225 
226  _check_size(_end_members);
227 
229 }
230 
232 template <typename... EndMembers>
234  EndMembers... _others) {
235  _add_end_member(_end_member);
236  _add_end_member(_others...);
237 }
238 
241 template <typename OutputIterator>
243  ParamComposition::AllowedOccupants _allowed_occs, OutputIterator result) {
244  ParamComposition param_comp(_allowed_occs);
245  param_comp.generate_prim_end_members();
246  param_comp.generate_composition_space();
247 
248  for (int i = 0; i < param_comp.allowed_list().size(); i++) {
249  const ParamComposition &curr = param_comp.allowed_list()[i];
250  Eigen::MatrixXd end_members(curr.origin().size(),
251  curr.spanning_end_members().size());
252  for (int j = 0; j < curr.spanning_end_members().size(); j++) {
253  end_members.col(j) = curr.spanning_end_members()[j];
254  }
255 
256  *result++ = CompositionConverter(param_comp.components().begin(),
257  param_comp.components().end(),
258  curr.origin(), end_members);
259  }
260 
261  return result;
262 }
263 
264 } // namespace CASM
265 
266 #endif
Convert between number of species per unit cell and parametric composition.
Eigen::VectorXd mol_composition(const Eigen::VectorXd &x) const
Convert parametric composition, 'x', to number of mol per prim, 'n'.
Eigen::MatrixXd m_end_members
Column vector matrix, rows == m_components.size(), cols == rank of parametric composition space.
std::string comp_formula(size_type i) const
Return formula for comp(i) in terms of comp_n(A), comp_n(B), ...
Eigen::VectorXd chem_pot(const Eigen::VectorXd param_chem_pot) const
Convert dG/dx to dG/dn.
std::string _n_formula(const Eigen::VectorXd &vec) const
Return formula for 'n'.
void _calc_conversion_matrices()
Calculate conversion matrices m_to_n and m_to_x.
Eigen::VectorXd m_origin
Vector, size == m_components.size(), specifying the num_mols_per_prim of each component at the origin...
Eigen::MatrixXd m_to_x
Conversion matrix: x = m_to_x*(n - origin)
Eigen::VectorXd end_member(size_type i) const
The mol composition of the parameteric composition axes end members.
std::string origin_formula() const
Return formula for origin.
std::string param_formula() const
Return formula for n->x.
void _check_size(const Eigen::MatrixXd &vec) const
Check that origin and end member vectors have same size as the number of components.
Eigen::VectorXd param_composition(const Eigen::VectorXd &n) const
Convert number of mol per prim, 'n' to parametric composition 'x'.
std::string comp_n_formula(size_type i) const
Return formula for comp_n(components()[i]) in terms of comp(a), comp(b), ...
static std::string comp_var(size_type i)
Composition variable names: "a", "b", ...
CompositionConverter()
Default constructor.
Eigen::VectorXd dparam_composition(const Eigen::VectorXd &dn) const
Convert change in number of atoms per prim, 'dn' to change in parametric composition 'dx'.
Eigen::VectorXd origin() const
The mol composition of the parameteric composition axes origin.
std::string param_chem_pot_formula(size_type i) const
Return formula for param_chem_pot(i) in terms of chem_pot(A), chem_pot(B), ...
size_type independent_compositions() const
The dimensionality of the composition space.
Eigen::VectorXd param_chem_pot(const Eigen::VectorXd chem_pot) const
Convert dG/dn to dG/dx.
std::string mol_formula() const
Return formula for x->n.
Eigen::VectorXd dmol_composition(const Eigen::VectorXd &dx) const
Convert change in parametric composition, 'dx', to change in number of mol per prim,...
void _add_end_member(Eigen::VectorXd _end_member)
Helps make variadic template constructor possible.
Eigen::MatrixXd dparam_dmol() const
Return the matrix Mij = dx_i/dn_j.
Eigen::MatrixXd dmol_dparam() const
Return the matrix Mij = dn_i/dx_j.
Eigen::MatrixXd m_to_n
Conversion matrix: n = origin + m_to_n*x.
std::vector< std::string > components() const
The order of components in mol composition vectors.
std::vector< std::string > m_components
List of all allowed components names in the prim, position in vector is reference for origin and end_...
std::string end_member_formula(size_type i) const
Return formula for end member.
const std::vector< ParamComposition > & allowed_list() const
std::vector< std::vector< std::string > > AllowedOccupants
const Eigen::VectorXd & origin() const
void generate_composition_space(bool verbose=false)
const std::vector< std::string > & components() const
Components are in order of appearance precedence in allowed_occs()
const std::vector< Eigen::VectorXd > & spanning_end_members() const
Main CASM namespace.
Definition: APICommand.hh:8
void display_comp_n(std::ostream &stream, const CompositionConverter &f, int indent=0)
Pretty-print comp_n in terms of comp.
Eigen::MatrixXd MatrixXd
Eigen::MatrixXd composition_space(const ParamComposition::AllowedOccupants &_allowed_occs, double tol=1e-14)
Return the composition space of a ParamComposition::AllowedOccupants.
void display_param_chem_pot(std::ostream &stream, const CompositionConverter &f, int indent=0)
Pretty-print param_chem_pot in terms of chem_pot.
void display_composition_axes(std::ostream &stream, const std::map< std::string, CompositionConverter > &map)
Pretty-print map of name/CompositionConverter pairs.
Eigen::MatrixXd end_members(const ParamComposition::AllowedOccupants &_allowed_occs)
Serialize CompositionConverter to JSON.
Eigen::MatrixXd null_composition_space(const ParamComposition::AllowedOccupants &_allowed_occs, double tol=1e-14)
Return the null composition space of a ParamComposition::AllowedOccupants.
void display_comp(std::ostream &stream, const CompositionConverter &f, int indent=0)
Pretty-print comp in terms of comp_n.
Eigen::VectorXd VectorXd
OutputIterator standard_composition_axes(ParamComposition::AllowedOccupants _allowed_occs, OutputIterator result)
Generate CompositionConverter specifying standard composition axes for a prim Structure.