CASM  1.1.0
A Clusters Approach to Statistical Mechanics
ConfigDoFValues.hh
Go to the documentation of this file.
1 #ifndef CASM_ConfigDoFValues
2 #define CASM_ConfigDoFValues
3 
6 namespace CASM {
7 
9  public:
10  ConfigDoFValues(DoF::BasicTraits const &_traits, Index _n_sublat,
11  Index _n_vol);
12 
13  std::string const &type_name() const;
14 
15  Index n_vol() const;
16 
17  Index n_sublat() const;
18 
19  private:
23 };
24 
39  public:
40  typedef Eigen::VectorXi ValueType;
41  typedef Eigen::VectorXi &Reference;
42  typedef Eigen::VectorXi const &ConstReference;
43 
44  typedef typename ValueType::Scalar SiteValueType;
45  typedef int &SiteReference;
46  typedef const int &ConstSiteReference;
47 
49  typedef typename ValueType::SegmentReturnType SublatReference;
50  typedef typename ValueType::ConstSegmentReturnType ConstSublatReference;
51 
52  LocalDiscreteConfigDoFValues(DoF::BasicTraits const &_traits, Index _n_sublat,
53  Index _n_vol,
54  std::vector<SymGroupRepID> const &_symrep_IDs);
55 
57  int &occ(Index i);
58 
60  void set_values(Eigen::Ref<ValueType const> const &_values);
61 
63  void setZero();
64 
67  Eigen::VectorXi const &values() const;
68 
71 
74 
77  std::vector<SymGroupRepID> const &symrep_IDs() const;
78 
79  private:
80  void _throw_if_invalid_size(Eigen::Ref<ValueType const> const &_values) const;
81 
83  std::vector<SymGroupRepID> m_symrep_IDs;
84 };
85 
107  public:
111 
113  typedef typename ValueType::ColXpr SiteReference;
114  typedef const typename ValueType::ConstColXpr ConstSiteReference;
115 
117  typedef typename Eigen::Block<ValueType> SublatReference;
118  typedef const typename Eigen::Block<const ValueType> ConstSublatReference;
119 
121  static Index matrix_dim(std::vector<DoFSetInfo> const &_info);
122 
124  Index _n_sublat, Index _n_vol,
125  std::vector<DoFSetInfo> const &_info);
126 
128  Index dim() const;
129 
131  void set_values(Eigen::Ref<const ValueType> const &_values);
132 
134  void setZero();
135 
137  Eigen::MatrixXd const &values() const;
138 
141  Eigen::Ref<const Eigen::MatrixXd> const &_standard_values);
142 
145 
148 
151 
154 
157 
159  std::vector<DoFSetInfo> const &info() const;
160 
161  private:
162  void _throw_if_invalid_size(Eigen::Ref<ValueType const> const &_values) const;
163 
165  std::vector<DoFSetInfo> m_info;
167 };
168 
169 // TODO: It might be confusing that ValueType, Reference, and ConstReference are
170 // in terms of Eigen::VectorXd, but dim(), from_standard_values(), and
171 // standard_values() are written using Eigen::MatrixXd. Since Eigen::VectorXd is
172 // a 1 column Eigen::MatrixXd it could be all written in terms of
173 // Eigen::VectorXd or Eigen::MatrixXd. Not sure if it is more useful to express
174 // in terms of Eigen::MatrixXd and match the LocalContinuousConfigDoFValues
175 // interface or use Eigen::VectorXd consistently to express that the value is 1
176 // dimensional.
177 
191  public:
195 
196  typedef typename ValueType::Scalar SiteValueType;
197  typedef int &SiteReference;
198  typedef const int &ConstSiteReference;
199 
201  Index _n_sublat, Index _n_vol,
202  DoFSetInfo const &_info);
203 
205  Index dim() const;
206 
208  void set_values(Eigen::Ref<const Eigen::MatrixXd> const &_values);
209 
211  void setZero();
212 
214  Eigen::VectorXd const &values() const;
215 
218  Eigen::Ref<const Eigen::MatrixXd> const &_standard_values);
219 
222 
224  DoFSetInfo const &info() const;
225 
226  private:
227  void _throw_if_invalid_size(Eigen::Ref<ValueType const> const &_values) const;
228 
231 };
232 
233 } // namespace CASM
234 
235 #endif
Specifies traits of (possibly) anisotropic crystal properties.
std::string const & type_name() const
ConfigDoFValues(DoF::BasicTraits const &_traits, Index _n_sublat, Index _n_vol)
Eigen::MatrixXd standard_values() const
Get global DoF values as standard DoF values.
GlobalContinuousConfigDoFValues(DoF::BasicTraits const &_traits, Index _n_sublat, Index _n_vol, DoFSetInfo const &_info)
void from_standard_values(Eigen::Ref< const Eigen::MatrixXd > const &_standard_values)
Set global DoF values from standard DoF values.
void _throw_if_invalid_size(Eigen::Ref< ValueType const > const &_values) const
Eigen::VectorXd const & values() const
Const access global DoF values.
Index dim() const
Global DoF vector representation dimension.
void setZero()
Set DoF values to zero.
DoFSetInfo const & info() const
DoFSetInfo provides the basis and symmetry representations for values
void set_values(Eigen::Ref< const Eigen::MatrixXd > const &_values)
Set global DoF values.
Eigen::VectorXd const & ConstReference
std::vector< DoFSetInfo > const & info() const
DoFSetInfo provides the basis and symmetry representations for values
SiteReference site_value(Index l)
Access site DoF value vector.
Eigen::MatrixXd const & ConstReference
Index dim() const
maximum DoF vector representation size (max of DoFSetInfo::dim())
const Eigen::Block< const ValueType > ConstSublatReference
static Index matrix_dim(std::vector< DoFSetInfo > const &_info)
local continuous DoF values matrix has #rows == max( DoFSetInfo::dim() )
void from_standard_values(Eigen::Ref< const Eigen::MatrixXd > const &_standard_values)
Set local DoF values from standard DoF values.
std::vector< DoFSetInfo > m_info
void _throw_if_invalid_size(Eigen::Ref< ValueType const > const &_values) const
Eigen::Block< ValueType > SublatReference
LocalContinuousConfigDoFValues(DoF::BasicTraits const &_traits, Index _n_sublat, Index _n_vol, std::vector< DoFSetInfo > const &_info)
const ValueType::ConstColXpr ConstSiteReference
void set_values(Eigen::Ref< const ValueType > const &_values)
Access site DoF values (prim DoF basis, matrix representing all sites)
SublatReference sublat(Index b)
Access matrix block of values for all sites on one sublattice.
Eigen::MatrixXd const & values() const
Const access DoF values (prim DoF basis, matrix representing all sites)
void setZero()
Set DoF values to zero.
Eigen::MatrixXd standard_values() const
Get local DoF values as standard DoF values.
std::vector< SymGroupRepID > const & symrep_IDs() const
ValueType::ConstSegmentReturnType ConstSublatReference
Eigen::VectorXi const & ConstReference
int & occ(Index i)
Reference occupation value on site i.
std::vector< SymGroupRepID > m_symrep_IDs
LocalDiscreteConfigDoFValues(DoF::BasicTraits const &_traits, Index _n_sublat, Index _n_vol, std::vector< SymGroupRepID > const &_symrep_IDs)
Eigen::VectorXi const & values() const
SublatReference sublat(Index b)
Access vector block of values for all sites on one sublattice.
void _throw_if_invalid_size(Eigen::Ref< ValueType const > const &_values) const
ValueType::SegmentReturnType SublatReference
void set_values(Eigen::Ref< ValueType const > const &_values)
Set occupation values (values are indices into Site::occupant_dof())
void setZero()
Set occupation values to zero.
Main CASM namespace.
Definition: APICommand.hh:8
Eigen::MatrixXd MatrixXd
std::string DoFKey
Definition: DoFDecl.hh:7
INDEX_TYPE Index
For long integer indexing:
Definition: definitions.hh:39
Eigen::VectorXd VectorXd