CASM
AClustersApproachtoStatisticalMechanics
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules
CompositionConverter.hh
Go to the documentation of this file.
1 #ifndef CASM_CompositionConverter_HH
2 #define CASM_CompositionConverter_HH
3 
4 #include <vector>
5 #include <map>
6 #include "casm/external/Eigen/Dense"
7 
10 
11 namespace CASM {
12 
13  class Structure;
14 
20 
21  public:
22 
23  typedef unsigned int size_type;
24 
27 
29  template<typename ComponentIterator, typename... EndMembers>
30  CompositionConverter(ComponentIterator begin, ComponentIterator end, Eigen::VectorXd _origin, EndMembers... _end_members);
31 
33  template<typename ComponentIterator>
34  CompositionConverter(ComponentIterator begin, ComponentIterator end, Eigen::VectorXd _origin, Eigen::MatrixXd _end_members);
35 
36 
38  size_type independent_compositions() const;
39 
41  static std::string comp_var(size_type i);
42 
44  std::vector<std::string> components() const;
45 
47  Eigen::VectorXd origin() const;
48 
50  Eigen::VectorXd end_member(size_type i) const;
51 
54 
57 
60 
63 
66 
69 
72 
75 
76 
78  std::string mol_formula() const;
79 
81  std::string param_formula() const;
82 
84  std::string origin_formula() const;
85 
87  std::string end_member_formula(size_type i) const;
88 
90  std::string comp_formula(size_type i) const;
91 
93  std::string comp_n_formula(size_type i) const;
94 
96  std::string param_chem_pot_formula(size_type i) const;
97 
98 
99  private:
100 
102  void _add_end_member(Eigen::VectorXd _end_member);
103 
105  template<typename... EndMembers>
106  void _add_end_member(Eigen::VectorXd _end_member, EndMembers... _others);
107 
109  void _check_size(const Eigen::VectorXd &vec) const;
110 
113 
115  std::string _n_formula(const Eigen::VectorXd &vec) const;
116 
117 
120  std::vector<std::string> m_components;
121 
125 
129 
133 
137 
138  };
139 
141  template<typename OutputIterator>
142  OutputIterator standard_composition_axes(const Structure &prim, OutputIterator result);
143 
145  void display_composition_axes(std::ostream &stream, const std::map<std::string, CompositionConverter> &map);
146 
148  void display_comp(std::ostream &stream, const CompositionConverter &f, int indent = 0);
149 
151  void display_comp_n(std::ostream &stream, const CompositionConverter &f, int indent = 0);
152 
154  void display_param_chem_pot(std::ostream &stream, const CompositionConverter &f, int indent = 0);
155 
158 
160  void from_json(CompositionConverter &f, const jsonParser &json);
161 
164 
166  Eigen::MatrixXd composition_space(const Structure &prim, double tol = 1e-14);
167 
169  Eigen::MatrixXd null_composition_space(const Structure &prim, double tol = 1e-14);
170 
171 
172  // ------ Definitions ---------------------------------------------
173 
183  template<typename ComponentIterator, typename... EndMembers>
185  ComponentIterator end,
186  Eigen::VectorXd _origin,
187  EndMembers... _end_members) :
188  m_components(begin, end),
189  m_origin(_origin) {
190 
191  _add_end_member(_end_members...);
192 
193  _check_size(_origin);
194 
196 
197  }
198 
208  template<typename ComponentIterator>
210  ComponentIterator end,
211  Eigen::VectorXd _origin,
212  Eigen::MatrixXd _end_members) :
213  m_components(begin, end),
214  m_origin(_origin),
215  m_end_members(_end_members) {
216 
217  _check_size(_origin);
218 
219  _check_size(_end_members.col(0));
220 
222 
223  }
224 
226  template<typename... EndMembers>
227  void CompositionConverter::_add_end_member(Eigen::VectorXd _end_member, EndMembers... _others) {
228  _add_end_member(_end_member);
229  _add_end_member(_others...);
230  }
231 
233  template<typename OutputIterator>
234  OutputIterator standard_composition_axes(const Structure &prim, OutputIterator result) {
235  ParamComposition param_comp(prim);
236  param_comp.generate_components();
237  param_comp.generate_sublattice_map();
238  param_comp.generate_prim_end_members();
239  param_comp.generate_composition_space();
240 
241  std::vector<std::string> components;
242  for(int i = 0; i < param_comp.get_components().size(); i++) {
243  components.push_back(param_comp.get_components()[i]);
244  }
245 
246  for(int i = 0; i < param_comp.get_allowed_list().size(); i++) {
247  const ParamComposition &curr = param_comp.get_allowed_list()[i];
249  for(int j = 0; j < curr.get_spanning_end_members().size(); j++) {
250  end_members.col(j) = curr.get_spanning_end_members()[j];
251  }
252 
253  *result++ = CompositionConverter(components.begin(), components.end(), curr.get_origin(), end_members);
254  }
255 
256  return result;
257  }
258 
259 
260 }
261 
262 #endif
263 
Eigen::MatrixXd MatrixXd
Eigen::VectorXd mol_composition(const Eigen::VectorXd &x) const
Convert parametric composition, 'x', to number of mol per prim, 'n'.
Eigen::MatrixXd dparam_dmol() const
Return the matrix Mij = dx_i/dn_j.
void display_comp_n(std::ostream &stream, const CompositionConverter &f, int indent=0)
Pretty-print comp_n in terms of comp.
void display_composition_axes(std::ostream &stream, const std::map< std::string, CompositionConverter > &map)
Pretty-print map of name/CompositionConverter pairs.
Eigen::MatrixXd composition_space(const Structure &prim, double tol=1e-14)
Return the composition space of a Structure.
void from_json(ClexDescription &desc, const jsonParser &json)
Index size() const
Definition: Array.hh:145
Eigen::MatrixXd dmol_dparam() const
Return the matrix Mij = dn_i/dx_j.
Eigen::VectorXd origin() const
The mol composition of the parameteric composition axes origin.
Eigen::VectorXd param_composition(const Eigen::VectorXd &n) const
Convert number of mol per prim, 'n' to parametric composition 'x'.
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), ...
Structure specifies the lattice and atomic basis of a crystal.
Definition: Structure.hh:29
jsonParser & to_json(const ClexDescription &desc, jsonParser &json)
std::vector< std::string > m_components
List of all allowed components names in the prim, position in vector is reference for origin and end_...
void _calc_conversion_matrices()
Calculate conversion matrices m_to_n and m_to_x.
void display_param_chem_pot(std::ostream &stream, const CompositionConverter &f, int indent=0)
Pretty-print param_chem_pot in terms of chem_pot.
const Array< ParamComposition > & get_allowed_list() const
std::string end_member_formula(size_type i) const
Return formula for end member.
Main CASM namespace.
Definition: complete.cpp:8
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, 'dn'.
Eigen::VectorXd chem_pot(const Eigen::VectorXd param_chem_pot) const
Convert dG/dx to dG/dn.
void _add_end_member(Eigen::VectorXd _end_member)
Helps make variadic template constructor possible.
std::vector< std::string > components() const
The order of components in mol composition vectors.
double tol
static std::string comp_var(size_type i)
Composition variable names: "a", "b", ...
Eigen::VectorXd param_chem_pot(const Eigen::VectorXd chem_pot) const
Convert dG/dn to dG/dx.
std::string comp_formula(size_type i) const
Return formula for comp(i) in terms of comp_n(A), comp_n(B), ...
Eigen::MatrixXd m_to_n
Conversion matrix: n = origin + m_to_n*x.
void _check_size(const Eigen::VectorXd &vec) const
Check that origin and end member vectors have same size as the number of components.
Eigen::MatrixXd m_end_members
Column vector matrix, rows == m_components.size(), cols == rank of parametric composition space...
Eigen::VectorXd VectorXd
const Array< Eigen::VectorXd > & get_spanning_end_members() const
size_type independent_compositions() const
The dimensionality of the composition space.
Eigen::VectorXd m_origin
Vector, size == m_components.size(), specifying the num_mols_per_prim of each component at the origin...
void generate_composition_space(bool verbose=false)
std::string origin_formula() const
Return formula for origin.
Eigen::MatrixXd m_to_x
Conversion matrix: x = m_to_x*(n - origin)
std::string comp_n_formula(size_type i) const
Return formula for comp_n(components()[i]) in terms of comp(a), comp(b), ...
OutputIterator standard_composition_axes(const Structure &prim, OutputIterator result)
Generate CompositionConverter specifying standard composition axes for a prim Structure.
Eigen::MatrixXd end_members(const Structure &prim)
Generate a column matrix containing all the possible molecular end members.
Convert between number of species per unit cell and parametric composition.
std::string param_formula() const
Return formula for n->x.
const Eigen::VectorXd & get_origin() const
void display_comp(std::ostream &stream, const CompositionConverter &f, int indent=0)
Pretty-print comp in terms of comp_n.
Eigen::VectorXd dparam_composition(const Eigen::VectorXd &dn) const
Convert change in number of atoms per prim, 'dn' to change in parametric composition 'dx'...
std::string _n_formula(const Eigen::VectorXd &vec) const
Return formula for 'n'.
Eigen::VectorXd end_member(size_type i) const
The mol composition of the parameteric composition axes end members.
CompositionConverter()
Default constructor.
Eigen::MatrixXd null_composition_space(const Structure &prim, double tol=1e-14)
Return the null composition space of a Structure.
const Array< std::string > & get_components() const
Components are ordered as in Structure::get_struc_molecule.