CASM
AClustersApproachtoStatisticalMechanics
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules
VaspIO.hh
Go to the documentation of this file.
1 #ifndef CASM_VaspIO
2 #define CASM_VaspIO
3 
7 #include "casm/clex/Supercell.hh"
9 
10 namespace CASM {
11 
13  namespace VaspIO {
14 
23  class SelectiveDynamics {
28 
29  public:
30 
33  m_option[0] = false;
34  m_option[1] = false;
35  m_option[2] = false;
36  }
37 
38  SelectiveDynamics(bool x, bool y, bool z) {
39  m_option[0] = x;
40  m_option[1] = y;
41  m_option[2] = z;
42  }
43 
44  bool &operator[](int i) {
45  return m_option[i];
46  }
47 
48  const bool &operator[](int i) const {
49  return m_option[i];
50  }
51 
52 
53  private:
54 
55  bool m_option[3];
56 
57  };
58 
63  inline std::ostream &operator<<(std::ostream &sout, const SelectiveDynamics &sel) {
64 
65  auto print = [&](bool val) {
66  if(val) {
67  sout << "T";
68  }
69  else {
70  sout << "F";
71  }
72  };
73 
74  print(sel[0]);
75  sout << " ";
76  print(sel[1]);
77  sout << " ";
78  print(sel[2]);
79  return sout;
80  };
81 
82 
83  namespace vaspio_impl {
84 
86 
87  public:
88 
101  m_title(""),
102  m_scale(1.0),
104  m_atom_names(true),
105  m_sel_dynamics(false),
106  m_append_atom_names(true),
107  m_ignore {"VA", "Va", "va"} {}
108 
118  PrintPOSCARBase(const Lattice &lat) :
119  PrintPOSCARBase() {
120  m_lat = lat;
121  }
122 
124  void set_title(std::string title) {
125  m_title = title;
126  }
127 
129  void set_scale(double s) {
130  m_scale = s;
131  }
132 
134  void set_direct() {
135  m_coord_mode = FRAC;
136  }
138  void set_frac() {
139  m_coord_mode = FRAC;
140  }
141 
143  void set_cart() {
144  m_coord_mode = CART;
145  }
146 
149  m_coord_mode = mode;
150  }
151 
154  m_sel_dynamics = false;
155  }
156 
159  m_sel_dynamics = true;
160  }
161 
164  m_atom_names = false;
165  }
166 
169  m_atom_names = true;
170  }
171 
174  m_append_atom_names = false;
175  }
176 
179  m_append_atom_names = false;
180  }
181 
183  std::vector<std::string> &ignore() {
184  return m_ignore;
185  }
186 
188  const std::vector<std::string> &ignore() const {
189  return m_ignore;
190  }
191 
192 
193  protected:
194 
196  template<typename TupleIterator>
197  void _print(std::ostream &sout, TupleIterator begin, TupleIterator end);
198 
199  private:
200 
201  std::string m_title;
202  double m_scale;
208 
210  std::vector<std::string> m_ignore;
211 
212  };
213 
214  }
215 
233 
234  public:
235 
236  typedef std::string AtomName;
237 
239  typedef std::tuple<AtomName, Coordinate, SelectiveDynamics> tuple_type;
240 
241  typedef std::vector<tuple_type>::iterator iterator;
242 
243  typedef std::vector<tuple_type>::const_iterator const_iterator;
244 
245 
247  PrintPOSCAR(const BasicStructure<Site> &struc);
248 
250  PrintPOSCAR(const Configuration &config);
251 
253  PrintPOSCAR(const Supercell &scel, const ConfigDoF &configdof);
254 
256  iterator begin() {
257  return m_atom_order.begin();
258  }
259 
261  iterator end() {
262  return m_atom_order.end();
263  }
264 
266  const_iterator cbegin() const {
267  return m_atom_order.cbegin();
268  }
269 
271  const_iterator cend() const {
272  return m_atom_order.cend();
273  }
274 
276  void sort();
277 
279  void print(std::ostream &sout);
280 
281 
282  private:
283 
285  std::vector<tuple_type> m_atom_order;
286 
287  };
288 
289 
290  // --- Definitions -------------------------------------------------------- //
291 
292  namespace vaspio_impl {
293 
295  template<typename TupleIterator>
296  void PrintPOSCARBase::_print(std::ostream &sout,
297  TupleIterator begin,
298  TupleIterator end) {
299 
300  int tprec = sout.precision();
301  std::ios::fmtflags tflags = sout.flags();
302  sout.precision(8);
303 
304  typedef std::tuple<std::string, Coordinate, SelectiveDynamics> tuple_type;
305 
306  // first filter out all atoms we are going to ignore, the remaining atoms get put in 'atom'
307  std::vector<tuple_type> atom;
308  for(auto it = begin; it != end; ++it) {
309 
310  // if Atom's name is not found in the ignore list, add it to 'atom'
311  if(m_ignore.cend() == std::find_if(m_ignore.cbegin(),
312  m_ignore.cend(),
313  [&](const std::string & name) {
314  if(std::get<0>(*it) == name) {
315  return true;
316  }
317  return false;
318  })) {
319  atom.push_back(*it);
320  }
321  }
322 
323  // print title, scale, and lattice
324  sout << m_title << "\n";
325  sout << std::fixed << std::setprecision(8) << m_scale << "\n";
326 
327  sout.flags(std::ios::showpoint | std::ios::fixed | std::ios::right);
328 
329  sout << ' ' << std::setw(16) << m_lat[0].transpose() << '\n';
330  sout << ' ' << std::setw(16) << m_lat[1].transpose() << '\n';
331  sout << ' ' << std::setw(16) << m_lat[2].transpose() << '\n';
332 
333  // if after filtering out ignored atoms none are left, return
334  if(atom.size() == 0) {
335  return;
336  }
337 
338  // count number each atom, and optionally print atom names line
339  std::vector<int> atom_count = {1};
340  auto it = atom.cbegin();
341  std::string curr_atom = std::get<0>(*it);
342  if(m_atom_names) {
343  sout << curr_atom << " ";
344  }
345  ++it;
346  for(; it != atom.cend(); ++it) {
347  if(std::get<0>(*it) != curr_atom) {
348  atom_count.push_back(1);
349  curr_atom = std::get<0>(*it);
350  if(m_atom_names) {
351  sout << curr_atom << " ";
352  }
353  }
354  else {
355  atom_count.back()++;
356  }
357  }
358  if(m_atom_names) {
359  sout << "\n";
360  }
361 
362  for(int i = 0; i < atom_count.size(); i++) {
363  sout << atom_count[i] << " ";
364  }
365  sout << "\n";
366 
367  // print 'Selective Dynamics' if using selective dynamics
368  if(m_sel_dynamics) {
369  sout << "Selective Dynamics\n";
370  }
371 
372  // print coord mode
373  sout << COORD_MODE::NAME(m_coord_mode) << "\n";
374 
375  // print all coordinates, and seletive dynamics settings, and atom names if applicable
376  for(auto it = atom.cbegin(); it != atom.cend(); ++it) {
377  std::get<1>(*it).print(sout, m_coord_mode);
378  if(m_sel_dynamics) {
379  sout << " " << std::get<2>(*it);
380  }
381  if(m_append_atom_names) {
382  sout << " " << std::get<0>(*it);
383  }
384  sout << "\n";
385  }
386  sout << "\n";
387 
388  sout.precision(tprec);
389  sout.flags(tflags);
390 
391  }
392  }
393 
411  vaspio_impl::PrintPOSCARBase(struc.lattice()) {
412 
413  set_title(struc.title);
414 
415  // create tuples collecting (Atom name, Coordinate, SelectiveDynamics) for each site
416  for(int i = 0; i < struc.basis.size(); ++i) {
417  m_atom_order.push_back(
418  tuple_type(
419  struc.basis[i].occ_name(),
420  struc.basis[i],
422  )
423  );
424  }
425  }
426 
442  inline PrintPOSCAR::PrintPOSCAR(const Configuration &config) :
443  vaspio_impl::PrintPOSCARBase(config.get_supercell().get_real_super_lattice()) {
444 
445  set_title(config.name());
446 
447  const Supercell &scel = config.get_supercell();
448 
449  // create tuples collecting (Atom name, Coordinate, SelectiveDynamics) for each site
450  for(int i = 0; i < config.size(); ++i) {
451  m_atom_order.push_back(
452  tuple_type(
453  config.get_mol(i).name,
454  scel.coord(i), // no displacement yet
456  )
457  );
458  }
459 
460  }
461 
477  inline PrintPOSCAR::PrintPOSCAR(const Supercell &scel, const ConfigDoF &configdof) :
478  vaspio_impl::PrintPOSCARBase(scel.get_real_super_lattice()) {
479 
480  // get occupant name for site i in configdof
481  auto occ_name = [&](int i) {
482  return scel.get_prim().basis[scel.get_b(i)].site_occupant()[configdof.occ(i)].name;
483  };
484 
485  // create tuples collecting (Atom name, Coordinate, SelectiveDynamics) for each site
486  for(int i = 0; i < configdof.size(); ++i) {
487  m_atom_order.push_back(
488  tuple_type(
489  occ_name(i),
490  scel.coord(i), // no displacement yet
492  )
493  );
494  }
495 
496  }
497 
498 
500  inline void PrintPOSCAR::sort() {
501  std::sort(m_atom_order.begin(),
502  m_atom_order.end(),
503  [ = ](const tuple_type & A, const tuple_type & B) {
504  return std::get<0>(A) < std::get<0>(B);
505  });
506  }
507 
509  inline void PrintPOSCAR::print(std::ostream &sout) {
510  _print(sout, m_atom_order.begin(), m_atom_order.end());
511  }
512 
513  }
514 
515 
516 }
517 
518 #endif
SelectiveDynamics()
Default Constructor sets to false for all directions.
Definition: VaspIO.hh:32
PrintPOSCARBase(const Lattice &lat)
Construct PrintPOSCAR object.
Definition: VaspIO.hh:118
std::vector< std::string > m_ignore
List of atom names which should not be printed (primarily for vacancies)
Definition: VaspIO.hh:210
Coordinate coord(const UnitCellCoord &bijk) const
Definition: Supercell.cc:45
int & occ(Index i)
Definition: ConfigDoF.hh:61
const bool & operator[](int i) const
Definition: VaspIO.hh:48
const_iterator cbegin() const
Iterate over tuples of (AtomName, Coordinate, SelectiveDynamics)
Definition: VaspIO.hh:266
std::vector< tuple_type >::iterator iterator
Definition: VaspIO.hh:241
Index size() const
Definition: Array.hh:145
const Molecule & get_mol(Index site_l) const
Molecule on site l.
const Structure & get_prim() const
Definition: Supercell.cc:74
void set_atom_names_off()
Do not print atom names line.
Definition: VaspIO.hh:163
void set_atom_names_on()
Print atom names line.
Definition: VaspIO.hh:168
PrintPOSCAR(const BasicStructure< Site > &struc)
Construct PrintPOSCAR object.
Definition: VaspIO.hh:410
std::tuple< AtomName, Coordinate, SelectiveDynamics > tuple_type
Atom name, Coordinate, SelectiveDynamics.
Definition: VaspIO.hh:239
Main CASM namespace.
Definition: complete.cpp:8
Represents a supercell of the primitive parent crystal structure.
Definition: Supercell.hh:37
std::vector< tuple_type >::const_iterator const_iterator
Definition: VaspIO.hh:243
Index size() const
*** ACCESSORS ***
Definition: ConfigDoF.hh:53
void set_selective_dynamics_off()
Set selective dynamics off.
Definition: VaspIO.hh:153
void _print(std::ostream &sout, TupleIterator begin, TupleIterator end)
Print POSCAR, provide a range of std::tuple ...
Definition: VaspIO.hh:296
void print(std::ostream &sout)
Print POSCAR to stream.
Definition: VaspIO.hh:509
std::ostream & operator<<(std::ostream &sout, const SelectiveDynamics &sel)
Write SelectiveDynamics options.
Definition: VaspIO.hh:63
std::string AtomName
Definition: VaspIO.hh:236
iterator begin()
Iterate over tuples of (AtomName, Coordinate, SelectiveDynamics)
Definition: VaspIO.hh:256
void set_selective_dynamics_on()
Set selective dynamics on.
Definition: VaspIO.hh:158
void set_cart()
Set coordinate mode to Cartesian.
Definition: VaspIO.hh:143
void sort()
Default sort is by atom name.
Definition: VaspIO.hh:500
std::vector< tuple_type > m_atom_order
(AtomName, Coordinate, SelectiveDynamics)
Definition: VaspIO.hh:285
A container class for the different degrees of freedom a Configuration might have.
Definition: ConfigDoF.hh:27
Array< CoordType > basis
Lattice vectors that specifies periodicity of the crystal.
PrintPOSCARBase()
Construct PrintPOSCAR object.
Definition: VaspIO.hh:100
iterator end()
Iterate over tuples of (AtomName, Coordinate, SelectiveDynamics)
Definition: VaspIO.hh:261
void set_direct()
Set coordinate mode to Direct (fractional)
Definition: VaspIO.hh:134
std::string title
User-specified name of this Structure.
void set_frac()
Set coordinate mode to fractional (Direct)
Definition: VaspIO.hh:138
Store SelectiveDynamics options.
Definition: VaspIO.hh:27
std::vector< std::string > & ignore()
Access vector of atom names which should not be printed, such as for vacancies.
Definition: VaspIO.hh:183
std::string name() const
SCELV_A_B_C_D_E_F/i.
void set_title(std::string title)
Set title.
Definition: VaspIO.hh:124
Index size() const
Returns number of sites, NOT the number of primitives that fit in here.
void set_coord_mode(COORD_TYPE mode)
Set coordinate mode.
Definition: VaspIO.hh:148
static std::string NAME()
get a string with the name of the active mode
const_iterator cend() const
Iterate over tuples of (AtomName, Coordinate, SelectiveDynamics)
Definition: VaspIO.hh:271
Print POSCAR with formating options.
Definition: VaspIO.hh:232
void set_scale(double s)
Set scaling factor.
Definition: VaspIO.hh:129
const std::vector< std::string > & ignore() const
const Access vector of atom names which should not be printed, such as for vacancies ...
Definition: VaspIO.hh:188
Supercell & get_supercell() const
Get the Supercell for this Configuration.
void set_append_atom_names_on()
Append atom name to end of each coordinate line.
Definition: VaspIO.hh:178
Index get_b(Index i) const
Definition: Supercell.hh:259
std::string name
Definition: Molecule.hh:112
void set_append_atom_names_off()
Do not append atom name to end of each coordinate line.
Definition: VaspIO.hh:173
SelectiveDynamics(bool x, bool y, bool z)
Definition: VaspIO.hh:38
A Configuration represents the values of all degrees of freedom in a Supercell.