CASM  1.1.0
A Clusters Approach to Statistical Mechanics
SimpleStructureTools.cc
Go to the documentation of this file.
2 
3 #include <string>
4 
7 #include "casm/clex/Supercell.hh"
12 
13 namespace CASM {
14 
15 namespace clex_SimpleStructureTools_impl {
16 
20 void _apply_dofs(xtal::SimpleStructure &_sstruc, ConfigDoF const &_config,
21  xtal::BasicStructure const &_reference,
22  std::vector<DoFKey> which_dofs);
23 
24 } // namespace clex_SimpleStructureTools_impl
25 
27  Supercell const &_scel, ConfigDoF const &_dof,
28  std::vector<DoFKey> const &_which_dofs) {
29  return make_simple_structure(_scel, _dof, MappedProperties(), _which_dofs,
30  false);
31 }
32 
34  Configuration const &_config, std::vector<DoFKey> const &_which_dofs,
35  bool relaxed) {
36  if (relaxed && is_calculated(_config)) {
37  return make_simple_structure(_config.supercell(), _config.configdof(),
38  _config.calc_properties(), _which_dofs, true);
39  }
40  // else
41  return make_simple_structure(_config.supercell(), _config.configdof(),
42  MappedProperties(), _which_dofs, false);
43 }
44 
46  Supercell const &_scel, ConfigDoF const &_dof,
47  MappedProperties const &_props, std::vector<DoFKey> const &_which_dofs,
48  bool relaxed) {
49  xtal::SimpleStructure result;
50 
51  result.mol_info.resize(_dof.size());
52  if (!relaxed) {
53  result.lat_column_mat = _scel.lattice().lat_column_mat();
54  for (Index b = 0, l = 0; b < _dof.n_sublat(); ++b) {
55  for (Index v = 0; v < _dof.n_vol(); ++v, ++l) {
56  result.mol_info.cart_coord(l) = _scel.coord(l).const_cart();
57  }
58  }
59  } else {
60  result.lat_column_mat = _props.global.at("latvec");
61  result.mol_info.coords = _props.site.at("coordinate");
62  }
63 
64  for (Index b = 0, l = 0; b < _dof.n_sublat(); ++b) {
65  for (Index v = 0; v < _dof.n_vol(); ++v, ++l) {
66  Molecule const &mol = _scel.prim().basis()[b].occupant_dof()[_dof.occ(l)];
67 
68  // Fill up the molecule's SpeciesAttributes
69  for (auto const &attr : mol.attributes()) {
70  // Has this attribute been encountered yet??
71  auto it = result.mol_info.properties.find(attr.first);
73  if (it == result.mol_info.properties.end()) {
74  // Iterator now points to initialized matrix
75  it = result.mol_info.properties
76  .emplace(attr.first,
77  Eigen::MatrixXd::Zero(attr.second.traits().dim(),
78  _dof.size()))
79  .first;
80  }
81  it->second.col(l) = attr.second.value();
82  }
83 
84  // Record name
85  result.mol_info.names[l] = mol.name();
86  }
87  }
88 
90  _which_dofs);
91  return result;
92 }
93 
94 //***************************************************************************
95 
96 std::vector<std::set<Index> > mol_site_compatibility(
97  xtal::SimpleStructure const &sstruc, Configuration const &_config) {
98  std::vector<std::set<Index> > result;
99  result.reserve(sstruc.mol_info.names.size());
100  for (std::string const &sp : sstruc.mol_info.names) {
101  result.push_back({});
102  for (Index l = 0; l < _config.size(); ++l) {
103  if (_config.mol(l).name() == sp) {
104  result.back().insert(l);
105  }
106  }
107  }
108  return result;
109 }
110 
111 std::vector<std::set<Index> > atom_site_compatibility(
112  xtal::SimpleStructure const &sstruc, Configuration const &_config) {
113  std::vector<std::set<Index> > result;
114  result.reserve(sstruc.atom_info.names.size());
115  for (std::string const &sp : sstruc.atom_info.names) {
116  result.push_back({});
117  for (Index l = 0; l < _config.size(); ++l) {
118  if (_config.mol(l).contains(sp)) {
119  result.back().insert(l);
120  }
121  }
122  }
123  return result;
124 }
125 
126 namespace clex_SimpleStructureTools_impl {
127 
137  public:
139  TransformDirective(std::string const &_name);
140 
142  std::string const &name() const { return m_name; }
143 
146  bool operator<(TransformDirective const &_other) const;
147 
150  void transform(ConfigDoF const &_config,
151  xtal::BasicStructure const &_reference,
152  xtal::SimpleStructure &_struc) const;
153 
154  private:
156  void _accumulate_before(std::set<std::string> const &_queue,
157  std::set<std::string> &_result) const;
158 
160  void _accumulate_after(std::set<std::string> const &_queue,
161  std::set<std::string> &_result) const;
162 
163  std::string m_name;
164  std::set<std::string> m_before;
165  std::set<std::string> m_after;
166 
168 };
169 
170 void _apply_dofs(xtal::SimpleStructure &_sstruc, ConfigDoF const &_config,
171  xtal::BasicStructure const &_reference,
172  std::vector<DoFKey> which_dofs) {
173  std::set<TransformDirective> tformers({TransformDirective("atomize")});
174  if (which_dofs.empty()) {
175  for (std::string const &dof : continuous_local_dof_types(_reference))
176  which_dofs.push_back(dof);
177  for (std::string const &dof : global_dof_types(_reference))
178  which_dofs.push_back(dof);
179  }
180 
181  for (DoFKey const &dof : which_dofs) {
182  if (dof != "none" && dof != "occ") tformers.insert(dof);
183  }
184 
185  for (TransformDirective const &tformer : tformers) {
186  tformer.transform(_config, _reference, _sstruc);
187  }
188 }
189 
190 TransformDirective::TransformDirective(std::string const &_name)
191  : m_name(_name), m_traits_ptr(nullptr) {
192  if (name() != "atomize") {
194  _accumulate_before({_name}, m_before);
195  _accumulate_after({_name}, m_after);
196  if (m_after.count("atomize") == 0) m_before.insert("atomize");
197  }
198 }
199 
201  if (m_before.count(_other.name()) || _other.m_after.count(name())) {
202  return false;
203  }
204  if (m_after.count(_other.name()) || _other.m_before.count(name())) {
205  return true;
206  }
207  return name() < _other.name();
208 }
209 
211  std::set<std::string> const &_queue, std::set<std::string> &_result) const {
212  for (std::string const &el : _queue) {
213  if (el != name()) _result.insert(el);
214  if (el != "atomize")
215  _accumulate_before(AnisoValTraits(el).must_apply_before(), _result);
216  }
217 }
218 
220  std::set<std::string> const &_queue, std::set<std::string> &_result) const {
221  for (std::string const &el : _queue) {
222  if (el != name()) _result.insert(el);
223  if (el != "atomize")
224  _accumulate_after(AnisoValTraits(el).must_apply_after(), _result);
225  }
226 }
227 
229  xtal::BasicStructure const &_reference,
230  xtal::SimpleStructure &_struc) const {
231  if (m_traits_ptr) {
232  if (m_traits_ptr->val_traits().global())
233  _struc.properties[m_traits_ptr->name()] =
235  else {
236  _struc.mol_info.properties[m_traits_ptr->name()] =
238  }
239  m_traits_ptr->apply_dof(_dof, _reference, _struc);
240  } else {
241  _atomize(_struc, _dof.occupation(), _reference);
242  }
243 }
244 } // namespace clex_SimpleStructureTools_impl
245 } // namespace CASM
Specifies traits of (possibly) anisotropic crystal properties.
bool global() const
Returns true if type is global.
Index size() const
Number of sites in the ConfigDoF.
Definition: ConfigDoF.cc:11
Index n_vol() const
Integer volume of ConfigDoF.
Definition: ConfigDoF.cc:14
Index n_sublat() const
Number of sublattices in ConfigDoF.
Definition: ConfigDoF.cc:17
GlobalContinuousConfigDoFValues const & global_dof(DoFKey const &_key) const
Definition: ConfigDoF.cc:69
Eigen::VectorXi const & occupation() const
Const reference occupation values.
Definition: ConfigDoF.cc:56
int & occ(Index i)
Reference occupation value on site i.
Definition: ConfigDoF.cc:34
LocalContinuousConfigDoFValues const & local_dof(DoFKey const &_key) const
Definition: ConfigDoF.cc:124
const ConfigDoF & configdof() const
const Access the DoF
const Molecule & mol(Index site_l) const
Molecule on site l.
const Supercell & supercell() const
Get the Supercell for this Configuration.
Index size() const
Returns number of sites, NOT the number of primitives that fit in here.
Collection of all the traits specific to a DoF type.
Definition: DoFTraits.hh:59
AnisoValTraits const & val_traits() const
Definition: DoFTraits.hh:67
virtual void apply_dof(ConfigDoF const &_dof, BasicStructure const &_reference, SimpleStructure &_struc) const
Transforms SimpleSructure.
Definition: DoFTraits.cc:99
std::string const & name() const
Definition: DoFTraits.hh:69
Eigen::MatrixXd standard_values() const
Get global DoF values as standard DoF values.
Eigen::MatrixXd standard_values() const
Get local DoF values as standard DoF values.
const std::vector< xtal::Site > & basis() const
Definition: Structure.hh:102
Represents a supercell of the primitive parent crystal structure.
Definition: Supercell.hh:51
Coordinate coord(Index linear_index) const
Return the linear index corresponding to integral coordinates.
Definition: Supercell.cc:199
const Lattice & lattice() const
The super lattice.
Definition: Supercell.cc:239
const Structure & prim() const
Definition: Supercell.cc:113
TransformDirective(std::string const &_name)
consturct from transformation or DoF type name
std::string const & name() const
Name of DoFType or transformation.
bool operator<(TransformDirective const &_other) const
Compare with _other TransformDirective. Returns true if this TransformDirective has precedence.
void transform(ConfigDoF const &_config, xtal::BasicStructure const &_reference, xtal::SimpleStructure &_struc) const
Applies transformation to _struc using information contained in _config.
void _accumulate_after(std::set< std::string > const &_queue, std::set< std::string > &_result) const
Build m_after object by recursively traversing DoF dependencies.
void _accumulate_before(std::set< std::string > const &_queue, std::set< std::string > &_result) const
Build m_before object by recursively traversing DoF dependencies.
BasicStructure specifies the lattice and atomic basis of a crystal.
const vector_type & const_cart() const
user override to force const Access the Cartesian coordinate vector
Definition: Coordinate.hh:90
const Eigen::Matrix3d & lat_column_mat() const
3x3 matrix with lattice vectors as its columne
Definition: Lattice.hh:110
Class representing a Molecule.
Definition: Molecule.hh:93
Representation of a crystal of molecular and/or atomic occupants, and any additional properties....
std::map< std::string, Eigen::MatrixXd > properties
std::string const & name() const
Designated name of Molecule (may be unrelated to constituent species)
Definition: Molecule.hh:117
bool contains(std::string const &atom_name) const
Returns true of molecule contains atom of specified name.
Definition: Molecule.cc:80
std::map< std::string, SpeciesAttribute > const & attributes() const
Returns dictionary of all constituent attributes of the Molecule Does not include attributes associat...
Definition: Molecule.hh:145
std::vector< DoFKey > continuous_local_dof_types(BasicStructure const &_struc)
std::vector< DoFKey > global_dof_types(BasicStructure const &_struc)
Traits const & traits(std::string const &dof_key)
Lookup DoFType::Traits in the global dictionary.
Definition: DoFTraits.cc:46
void _apply_dofs(xtal::SimpleStructure &_sstruc, ConfigDoF const &_config, xtal::BasicStructure const &_reference, std::vector< DoFKey > which_dofs)
Imposes DoF values from ConfigDoF _config onto *this, using using any necessary information contained...
void _atomize(SimpleStructure &_sstruc, Eigen::Ref< const Eigen::VectorXi > const &_mol_occ, BasicStructure const &_reference)
use information in _reference to initialize atom_info from mol_info
Main CASM namespace.
Definition: APICommand.hh:8
std::vector< std::set< Index > > atom_site_compatibility(xtal::SimpleStructure const &sstruc, Configuration const &_config)
Determine which sites of a Configuration can host each atom of a SimpleStructure result[i] is set of ...
std::vector< std::set< Index > > mol_site_compatibility(xtal::SimpleStructure const &sstruc, Configuration const &_config)
Determine which sites of a Configuration can host each molecule of a SimpleStructure result[i] is set...
xtal::SimpleStructure make_simple_structure(Supercell const &_scel, ConfigDoF const &_dof, std::vector< DoFKey > const &_which_dofs={})
Construct from ConfigDoF _dof belonging to provided Supercell _scel.
bool is_calculated(const ConfigType &config, std::string calctype="")
Return true if all required properties have been been calculated for the configuration.
Definition: Calculable.cc:137
std::string DoFKey
Definition: DoFDecl.hh:7
INDEX_TYPE Index
For long integer indexing:
Definition: definitions.hh:39
std::map< std::string, Eigen::MatrixXd > site
std::map< std::string, Eigen::MatrixXd > global
Eigen::MatrixXd coords
(3 x names.size()) matrix of coordinates. coords.col(i) is Cartesian coordinate of site 'i'
std::vector< std::string > names
names[i] is name of species that occupies sites 'i'
Coord cart_coord(Index i)
Access coordinate of site 'i'.
std::map< std::string, Eigen::MatrixXd > properties
map of [property name, (m x names.size()) matrix] for all numerical site properties properties are as...
void resize(Index N)
Resize to hold N sites. All coordinates set to zero, all occupants set to "Va".