CASM  1.1.0
A Clusters Approach to Statistical Mechanics
super.cc
Go to the documentation of this file.
1 #include <boost/filesystem/fstream.hpp>
2 #include <fstream>
3 
26 #include "casm/external/Eigen/Core"
27 #include "casm/external/Eigen/src/Core/Matrix.h"
29 
30 namespace CASM {
31 
32 namespace Completer {
34 
35 const std::vector<fs::path> &SuperOption::transf_mat_paths() const {
36  return m_transf_mat_paths;
37 }
38 
39 const fs::path &SuperOption::struct_path() const { return m_struct_path; }
40 
41 const std::string &SuperOption::unit_scel_str() const {
42  return m_unit_scel_str;
43 }
44 
46 
47 double SuperOption::tolerance() const { return m_tolerance; }
48 
55 
56  m_desc.add_options()("transf-mat",
57  po::value<std::vector<fs::path> >(&m_transf_mat_paths)
58  ->multitoken()
59  ->value_name(ArgHandler::path()),
60  "1 or more files containing a 3x3 transformation matrix "
61  "used to create a supercell.")
62 
63  ("get-transf-mat", "If it exists, find the transformation matrix.")
64 
65  ("structure",
66  po::value<fs::path>(&m_struct_path)->value_name(ArgHandler::path()),
67  "File with structure (POSCAR type) to use.")
68 
69  ("unitcell",
70  po::value<std::string>(&m_unit_scel_str)
71  ->value_name(ArgHandler::supercell()),
72  "Name of supercell to use as unit cell. For ex. "
73  "'SCEL2_2_1_1_0_0_0'.")
74 
75  ("duper",
76  "Construct the superduperlattice, the minimum supercell of "
77  "all input supercells "
78  "and configurations.")
79 
80  ("fixed-orientation",
81  "When constructing the superduperlattice, do not "
82  "consider other symmetrically "
83  "equivalent orientations.")
84 
85  ("min-volume", po::value<Index>(&m_min_vol),
86  "Transforms the transformation matrix, T -> T', "
87  "where T' = T*M, such that "
88  "(T').determininant() >= V. This has the effect "
89  "that a supercell has a "
90  "particular volume.")
91 
92  ("fixed-shape",
93  "Used with --min-volume to enforce that T' = "
94  "T*m*I, where I is the identity "
95  "matrix, and m is a scalar. This has the effect "
96  "of preserving the shape "
97  "of the resulting supercell, but increasing the "
98  "volume.")
99 
100  ("verbose",
101  "When used with --duper, show how the input "
102  "lattices are transformed "
103  "to tile the superduperlattice.")
104 
105  ("add-canonical,a",
106  "Will add the generated super "
107  "configuration in it's "
108  "canonical form in the equivalent "
109  "niggli supercell.")
110 
111  ("vasp5",
112  "Print using VASP5 style (include "
113  "atom name line)")
114 
115  ("tol",
116  po::value<double>(&m_tolerance)
117  ->default_value(CASM::TOL),
118  "Tolerance used for checking "
119  "symmetry");
120 
121  return;
122 }
123 } // namespace Completer
124 
125 // ///////////////////////////////////////
126 // 'super' function for casm
127 // (add an 'if-else' statement in casm.cpp to call this)
128 
129 int super_command(const CommandArgs &args) {
130  // casm enum [—supercell min max] [—config supercell ] [—hopconfigs
131  // hop.background]
132  //- enumerate supercells and configs and hop local configurations
133 
134  std::vector<std::string> scelname, configname;
135  std::string unitscelname;
136  fs::path structfile;
137  std::vector<fs::path> tmatfile, abs_tmatfile, config_path;
138  fs::path abs_structfile;
139  Index min_vol;
140  double tol;
141  COORD_TYPE coordtype;
142  po::variables_map vm;
143 
145  Completer::SuperOption super_opt;
146 
147  try {
148  po::store(
149  po::parse_command_line(args.argc(), args.argv(), super_opt.desc()),
150  vm); // can throw
151 
152  if (!vm.count("help") && !vm.count("desc")) {
153  if (!vm.count("duper")) {
154  if (vm.count("transf-mat") + vm.count("get-transf-mat") != 1) {
155  err_log() << "Error in 'casm super'. Only one of --transf-mat or "
156  "--get-transf-mat may be chosen."
157  << std::endl;
158  return ERR_INVALID_ARG;
159  }
160  if (configname.size() > 1 || scelname.size() > 1 ||
161  tmatfile.size() > 1) {
162  err_log() << "ERROR: more than one --confignames, --scelnames, or "
163  "--transf-mat argument "
164  "is only allowed for option --duper"
165  << std::endl;
166  return ERR_INVALID_ARG;
167  }
168  if (config_path.size() > 0) {
169  err_log() << "ERROR: the --configs option is only allowed with "
170  "option --duper"
171  << std::endl;
172  return ERR_INVALID_ARG;
173  }
174  }
175  }
176 
179  if (vm.count("help")) {
180  log() << "\n";
181  log() << super_opt.desc() << std::endl;
182 
183  return 0;
184  }
185 
186  if (vm.count("desc")) {
187  log() << "\n";
188  log() << super_opt.desc() << std::endl;
189  log() << "DESCRIPTION" << std::endl;
190  log() << " "
191  " \n"
192  << " casm super --transf-mat T "
193  " \n"
194  << " - Print super lattice of the PRIM lattice "
195  " \n"
196  << " "
197  " \n"
198  << " casm super --structure POSCAR --transf-mat T "
199  " \n"
200  << " - Print superstructure of a POSCAR "
201  " \n"
202  << " "
203  " \n"
204  << " casm super --confignames configname --transf-mat T "
205  " \n"
206  << " - Print superstructure of a configuration "
207  " \n"
208  << " "
209  " \n"
210  << " casm super --structure POSCAR --unitcell scelname "
211  "--get-transf-mat \n"
212  << " - Check if POSCAR lattice is a supercell of unit cell "
213  "lattice and \n"
214  << " if so print the transformation matrix "
215  " \n"
216  << " - Uses primitive cell for unitcell if none given "
217  " \n"
218  << " "
219  " \n"
220  << " casm super --scelnames scelname --unitcell scelname "
221  "--get-transf-mat\n"
222  << " - Check if configuration lattice is a supercell of unit cell "
223  "lattice.\n"
224  << " and print the transformation matrix "
225  " \n"
226  << " - Uses primitive cell for unitcell if none given "
227  " \n\n"
228  <<
229 
230  " casm super --duper --scelnames scel1 [scel2 ...] --confignames "
231  "con1 [con2 ...]\n"
232  " --configs [mylist ...] --transf-mat M1 [M2 ...] "
233  " \n"
234  << " - Makes the superduperlattice of the lattices of all inputs "
235  " \n"
236  << " - Using '--configs' with no arguments is equivalent to "
237  "'--configs MASTER',\n"
238  << " which uses the master config list "
239  " \n"
240  << " - Default applies prim point group ops to try to find "
241  "minimum volume\n"
242  << " superduperlattice, disable with '--fixed-orientation' "
243  " \n\n";
244 
245  return 0;
246  }
247 
248  po::notify(vm); // throws on error, so do after help in case
249  // there are any problems
250 
251  scelname = super_opt.supercell_strs();
252  configname = super_opt.config_strs();
253  unitscelname = super_opt.unit_scel_str();
254  structfile = super_opt.struct_path();
255  tmatfile = super_opt.transf_mat_paths();
256  config_path = super_opt.selection_paths();
257  min_vol = super_opt.min_vol();
258  tol = super_opt.tolerance();
259  coordtype = super_opt.coordtype_enum();
260  } catch (po::error &e) {
261  err_log() << "ERROR: " << e.what() << std::endl << std::endl;
262  err_log() << super_opt.desc() << std::endl;
263  return 1;
264  } catch (std::exception &e) {
265  err_log() << "Unhandled Exception reached the top of main: " << e.what()
266  << ", application will now exit" << std::endl;
267  return 1;
268  }
269 
270  xtal::COORD_MODE C(coordtype);
271 
272  // lambda for printing
273  auto print = [&](const BasicStructure &struc) {
274  VaspIO::PrintPOSCAR printer(make_simple_structure(struc), struc.title());
275 
276  if (vm.count("vasp5")) {
277  printer.set_atom_names_on();
278  } else {
279  printer.set_atom_names_off();
280  }
281  printer.set_coord_mode(coordtype);
282  printer.print(log());
283  };
284 
285  // -- no casm project necessary for super cell of a POSCAR -------
286 
287  // want absolute paths
288  for (auto &&file : tmatfile) {
289  abs_tmatfile.push_back(fs::absolute(file));
290  }
291  abs_structfile = fs::absolute(structfile);
292 
293  if (vm.count("structure") && vm.count("transf-mat")) {
294  if (!fs::exists(abs_structfile)) {
295  log() << "ERROR: " << abs_tmatfile[0] << " not found." << std::endl;
296  return 1;
297  }
298  std::ifstream abs_structfile_stream(abs_structfile.string());
299  BasicStructure unitcell =
300  BasicStructure::from_poscar_stream(abs_structfile_stream);
301 
302  // -- read transf matrix ---
303  if (!fs::exists(abs_tmatfile[0])) {
304  log() << "ERROR: " << abs_tmatfile[0] << " not found." << std::endl;
305  return 1;
306  }
307  Eigen::Matrix3i Tm;
308  fs::ifstream file(abs_tmatfile[0]);
309  file >> Tm;
310  file.close();
311 
312  BasicStructure super = xtal::make_superstructure(unitcell, Tm);
313  super.set_title(std::string("Supercell of ") + unitcell.title());
314 
315  print(super);
316 
317  return 0;
318  }
319 
320  const fs::path &root = args.root;
321  if (root.empty()) {
322  err_log().error("No casm project found");
323  err_log() << std::endl;
324  return ERR_NO_PROJ;
325  }
326 
327  // If 'args.primclex', use that, else construct PrimClex in 'uniq_primclex'
328  // Then whichever exists, store reference in 'primclex'
329  std::unique_ptr<PrimClex> uniq_primclex;
330  PrimClex &primclex = make_primclex_if_not(args, uniq_primclex);
331 
332  if (vm.count("duper")) {
333  // collect all the Lattice to make the superduperlattice of
334  std::map<std::string, Lattice> lat;
335  std::map<std::string, Lattice> config_lat;
336 
337  // collect lattices by constructing from transformation matrices
338  if (vm.count("transf-mat")) {
339  for (auto it = abs_tmatfile.begin(); it != abs_tmatfile.end(); ++it) {
340  Eigen::Matrix<int, 3, 3, Eigen::RowMajor> T;
341  fs::ifstream file(*it);
342  for (int i = 0; i < 9; i++) {
343  file >> T.data()[i];
344  }
345  file.close();
346  lat[it->string()] = make_superlattice(primclex.prim().lattice(), T);
347  }
348  }
349 
350  // collect supercells from --scelnames
351  if (vm.count("scelnames")) {
352  for (auto it = scelname.begin(); it != scelname.end(); ++it) {
353  lat[*it] = primclex.db<Supercell>().find(*it)->lattice();
354  }
355  }
356 
357  // collect configs from --confignames
358  if (vm.count("confignames")) {
359  for (auto it = configname.begin(); it != configname.end(); ++it) {
360  config_lat[*it] = lat[*it] =
361  primclex.db<Configuration>().find(*it)->supercell().lattice();
362  }
363  }
364 
365  // collect configs from lists via --configs
366  if (vm.count("configs")) {
367  // MASTER config list if '--configs' only
368  if (config_path.size() == 0) {
369  config_path.push_back("MASTER");
370  }
371 
372  for (const auto &p : config_path) {
374  for (const auto &config : selection.selected()) {
375  config_lat[config.name()] = lat[config.name()] =
376  config.supercell().lattice();
377  }
378  }
379  }
380 
381  std::vector<Lattice> lat_only;
382  for (auto it = lat.begin(); it != lat.end(); ++it) {
383  lat_only.push_back(it->second);
384  }
385 
386  // create superduperlattice
387  auto begin = primclex.prim().point_group().begin();
388  auto end = primclex.prim().point_group().end();
389  if (vm.count("fixed-orientation")) {
390  end = begin;
391  }
393  lat_only.begin(), lat_only.end(), begin, end);
394 
396  if (vm.count("min-volume")) {
397  log() << " Enforcing minimum volume: " << min_vol;
398  if (vm.count("fixed-shape")) {
399  log() << " (with fixed shape)";
400  }
401  log() << "\n\n";
402 
403  auto prim_lat = primclex.prim().lattice();
404  const SymGroup &pg = primclex.prim().point_group();
405  auto T = is_superlattice(superduper, prim_lat, TOL).second;
406 
407  log() << " superduperlattice lattice: \n"
408  << superduper.lat_column_mat() << "\n\n";
409 
410  log() << " Initial transformation matrix:\n"
411  << iround(T)
412  << "\n (volume = " << iround(T).cast<double>().determinant()
413  << ")\n\n";
414 
415  auto M = enforce_min_volume(pg.begin(), pg.end(), prim_lat, iround(T),
416  min_vol, vm.count("fixed-shape"));
417 
418  superduper = xtal::canonical::equivalent(make_superlattice(superduper, M),
419  pg, TOL);
420 
421  auto S = is_superlattice(superduper, prim_lat, TOL).second;
422 
423  log() << " superduperlattice lattice: \n"
424  << superduper.lat_column_mat() << "\n\n";
425 
426  log() << " Transformation matrix, after enforcing mininum volume:\n"
427  << iround(S)
428  << "\n (volume = " << iround(S).cast<double>().determinant()
429  << ")\n\n";
430  }
431 
432  const Supercell &superduper_scel =
433  *Supercell(&primclex, superduper).insert().first;
434  superduper = superduper_scel.lattice();
435 
436  log() << "--- Lattices as column vector matrices ---\n\n";
437 
438  log() << " superduperlattice: " << superduper_scel.name() << "\n\n";
439 
440  log() << " superduperlattice lattice: \n"
441  << superduper.lat_column_mat() << "\n\n";
442 
443  log() << " Transformation matrix, relative the primitive cell:\n";
444  log() << iround(is_superlattice(superduper, primclex.prim().lattice(), TOL)
445  .second)
446  << "\n\n";
447 
448  if (vm.count("verbose")) {
449  log() << "Transformation matrices: \n";
450  for (auto it = lat.begin(); it != lat.end(); ++it) {
451  log() << "--- \n";
452  log() << " Unit: " << it->first << ":\n"
453  << it->second.lat_column_mat() << "\n\n";
454 
455  auto res = xtal::is_equivalent_superlattice(superduper, it->second,
456  begin, end, TOL);
457  log() << " Superduper = (op*unit) * T\n\nop:\n";
458  log() << res.first->matrix() << "\n\n";
459  log() << " T:\n";
460  log() << iround(res.second) << "\n\n";
461  }
462  log() << "--- \n";
463  }
464 
465  log() << "Write supercell database..." << std::endl;
467  log() << " DONE" << std::endl << std::endl;
468 
469  if (vm.count("add-canonical")) {
470  log() << "Add super configurations:\n";
471  for (auto it = config_lat.begin(); it != config_lat.end(); ++it) {
472  FillSupercell f{superduper_scel};
473  auto insert_res =
474  f(*primclex.db<Configuration>().find(it->first)).insert();
475  if (insert_res.insert_canonical) {
476  log() << " " << it->first << " -> "
477  << insert_res.canonical_it->name() << "\n";
478  }
479  }
480  log() << "\n";
481  log() << "Writing configuration database..." << std::endl;
483  log() << " DONE" << std::endl;
484  }
485 
486  return 0;
487 
488  } else if (vm.count("transf-mat")) {
489  Eigen::Matrix3i T;
490  if (!fs::exists(abs_tmatfile[0])) {
491  log() << "ERROR: " << abs_tmatfile[0] << " not found." << std::endl;
492  return 1;
493  }
494  fs::ifstream file(abs_tmatfile[0]);
495  file >> T;
496  file.close();
497 
498  log() << "Read transformation matrix, T: \n" << T << "\n\n";
499 
501  if (vm.count("min-volume")) {
502  if (!vm.count("fixed-shape")) {
503  log() << " Enforcing minimum volume: \n";
504  log() << " Finding T' = T*M, such that (T').determinant() >= "
505  << min_vol;
506  } else {
507  log() << " Enforcing minimum volume (with fixed shape): \n";
508  log() << " Finding T' = T*m*I, such that (T').determinant() >= "
509  << min_vol;
510  }
511  log() << "\n\n";
512 
513  auto prim_lat = primclex.prim().lattice();
514  const SymGroup &pg = primclex.prim().point_group();
515 
516  log() << " Initial transformation matrix:\n"
517  << T << "\n (volume = " << T.cast<double>().determinant()
518  << ")\n\n";
519 
520  auto M =
521  enforce_min_volume(pg.begin(), pg.end(), primclex.prim().lattice(), T,
522  min_vol, vm.count("fixed-shape"));
523 
524  Eigen::Matrix3i super_lat_matrix = T * M;
526  make_superlattice(prim_lat, super_lat_matrix), pg, TOL);
527  T = iround(is_superlattice(niggli_lat, prim_lat, TOL).second);
528 
529  log() << " Transformation matrix, after enforcing mininum volume:\n"
530  << T << "\n (volume = " << iround(T).cast<double>().determinant()
531  << ")\n\n";
532  }
533 
534  // super lattice
535  if (vm.count("scelnames")) {
536  const Supercell &scel = *primclex.db<Supercell>().find(scelname[0]);
537 
538  log() << " Unit cell: " << scelname[0] << "\n\n";
539 
540  log() << " Unit cell lattice: \n"
541  << scel.lattice().lat_column_mat() << "\n\n";
542 
543  Lattice super_lat = make_superlattice(scel.lattice(), T);
544  const Supercell &super_scel =
545  *Supercell(&primclex, super_lat).insert().first;
546 
547  log() << " Add supercell: " << super_scel.name() << "\n\n";
548 
549  log() << " Supercell lattice: \n"
550  << super_scel.lattice().lat_column_mat() << "\n\n";
551 
552  log() << " Transformation matrix: \n"
553  << super_scel.transf_mat() << "\n\n";
554 
555  log() << "Write supercell database..." << std::endl;
557  log() << " DONE" << std::endl << std::endl;
558 
559  }
560  // super structure
561  else if (vm.count("confignames")) {
562  std::stringstream ss;
563  const Configuration &con =
565 
566  VaspIO::PrintPOSCAR p(make_simple_structure(con), con.name());
567  p.sort();
568  p.print(ss);
569 
570  std::istringstream iss(ss.str());
572 
573  log() << "Unit structure:";
574  log() << "\n------\n";
575  print(unit);
576  log() << "\n------\n";
577  log() << "\n\n";
578 
579  // TODO: Why is there so much code copied and pasted in all the casm
580  // commands?
581  BasicStructure super = xtal::make_superstructure(unit, T);
582  super.set_title(std::string("Supercell of ") + con.name());
583 
584  log() << "Super structure:";
585  log() << "\n------\n";
586  print(super);
587  log() << "\n------\n";
588 
589  if (vm.count("add-canonical")) {
590  double tol = TOL;
591 
592  ConfigMapper configmapper(primclex, ConfigMapping::Settings(), tol);
593 
594  auto map_res =
595  configmapper.import_structure(make_simple_structure(super));
596 
597  if (map_res.success()) {
598  auto insert_res = (map_res.maps.begin()->second).config.insert();
599  Configuration imported_config = *insert_res.canonical_it;
600 
601  if (insert_res.insert_canonical) {
602  log() << " The configuration was imported successfully as "
603  << imported_config.name() << std::endl
604  << std::endl;
605 
606  } else {
607  log() << " The configuration was mapped onto pre-existing "
608  "equivalent structure "
609  << imported_config.name() << std::endl
610  << std::endl;
611  }
612  jsonParser json_src;
613  json_src["supercell_of"] = configname[0];
614  imported_config.push_back_source(json_src);
615  primclex.db<Configuration>().update(imported_config);
616  } else {
617  log() << " Failure to map configuration." << std::endl << std::endl;
618  }
619 
620  // Update directories
621  log() << "Write supercell database..." << std::endl;
623  log() << " DONE" << std::endl << std::endl;
624 
625  log() << "Writing configuration database..." << std::endl;
627  log() << " DONE" << std::endl;
628  }
629 
630  return 0;
631 
632  }
633  // super lattice of prim lattice
634  else {
635  BasicStructure unit = primclex.prim();
636  SymGroup pg = Structure(unit).point_group();
637 
639  Eigen::Matrix3d S = U * T.cast<double>();
640  Eigen::Matrix3i H_canon;
641  Eigen::Matrix3d op_canon;
642 
643  Eigen::Matrix3d S_niggli =
645  Eigen::Matrix3i T_niggli = iround(U.inverse() * S_niggli);
646  Eigen::Matrix3i H_niggli = hermite_normal_form(T_niggli).first;
647 
648  std::stringstream s_name;
649  s_name << "SCEL" << H_niggli(0, 0) * H_niggli(1, 1) * H_niggli(2, 2)
650  << "_" << H_niggli(0, 0) << "_" << H_niggli(1, 1) << "_"
651  << H_niggli(2, 2) << "_" << H_niggli(1, 2) << "_" << H_niggli(0, 2)
652  << "_" << H_niggli(0, 1);
653 
654  // S = U*T;
655  // S_canon = op_canon*S = U*T', where H_canon*V = U.inv*op_canon*U*T = T'
656 
657  log() << "--- Lattices as column vector matrices ---\n\n";
658 
659  log() << "Prim lattice, U:\n" << U << "\n\n";
660 
661  log() << "Super lattice, S = U*T:\n" << S << "\n\n";
662 
663  log() << "This is equivalent to '" << s_name.str()
664  << "', the equivalent super lattice \n"
665  << "in the standard orientation niggli cell, S_niggli:\n"
666  << S_niggli << "\n\n";
667 
668  log() << "The transformation matrix (S_niggli = U*T) for '"
669  << s_name.str() << "' is:\n"
670  << T_niggli << "\n\n";
671 
672  log() << "--- Lattices as row vector matrices ---\n\n";
673 
674  log() << "Prim lattice:\n" << U.transpose() << "\n\n";
675 
676  log() << "Super lattice:\n" << S.transpose() << "\n\n";
677 
678  log() << "This is equivalent to '" << s_name.str()
679  << "', the equivalent super lattice \n"
680  << "in the standard orientation niggli cell:\n"
681  << S_niggli.transpose() << "\n\n";
682 
683  return 0;
684  }
685  }
686 
687  if (vm.count("get-transf-mat")) {
688  Lattice unit_lat = primclex.prim().lattice();
689 
690  if (vm.count("unitcell")) {
691  unit_lat = primclex.db<Supercell>().find(unitscelname)->lattice();
692  }
693 
694  Lattice super_lat;
695 
696  if (vm.count("structure")) {
697  std::ifstream abs_structfile_stream(abs_structfile.string());
698  super_lat =
699  BasicStructure::from_poscar_stream(abs_structfile_stream).lattice();
700  } else if (vm.count("scelnames")) {
701  super_lat = primclex.db<Supercell>().find(unitscelname)->lattice();
702  } else {
703  log() << "Error in 'casm super --get-transf-mat'. No --structure or "
704  "--scelnames given."
705  << std::endl
706  << std::endl;
707  return 1;
708  }
709 
710  log() << "--- Lattices as column vector matrices ---\n\n";
711 
712  log() << "Unit lattice, U:\n" << unit_lat.lat_column_mat() << "\n\n";
713 
714  log() << "Super lattice, S:\n" << super_lat.lat_column_mat() << "\n\n";
715 
716  // see if super_lat is a supercell of unitlat
717  // S == U*T
718  Eigen::Matrix3d T;
719  bool is_superlat;
720  std::tie(is_superlat, T) =
721  xtal::is_superlattice(super_lat, unit_lat, super_lat.tol());
722  if (is_superlat) {
723  log() << "The super lattice is a supercell of the unit lattice.\n\n";
724 
725  log() << "The transformation matrix, T, where S = U*T, is: \n"
726  << iround(T) << "\n\n";
727  } else {
728  log() << "The super lattice is NOT a supercell of the unit lattice.\n\n";
729 
730  log() << "The transformation matrix, T, where S = U*T, is: \n"
731  << T << "\n\n";
732  }
733 
734  return 0;
735  }
736 
737  return 0;
738 };
739 
740 } // namespace CASM
#define ERR_NO_PROJ
Definition: errors.hh:13
#define ERR_INVALID_ARG
Definition: errors.hh:7
int argc() const
Definition: CLIParse.hh:20
char ** argv() const
Definition: CLIParse.hh:22
static std::string path()
Get value_type string for path completion.
Definition: Handlers.cc:56
static std::string supercell()
Get value_type string for supercell completion.
Definition: Handlers.cc:60
void add_confignames_suboption()
Add a –confignames suboption.
Definition: Handlers.cc:691
const po::options_description & desc()
Get the program options, filled with the initialized values.
Definition: Handlers.cc:321
void add_help_suboption()
Add a plain –help and –desc suboptions.
Definition: Handlers.cc:561
po::options_description m_desc
Definition: Handlers.hh:260
void add_coordtype_suboption()
Add a –coord suboption to specify FRAC or CART.
Definition: Handlers.cc:725
void add_scelnames_suboption()
Add a –scelnames suboption.
Definition: Handlers.cc:665
void add_configlists_nodefault_suboption()
Add –configs suboption (no default)
Definition: Handlers.cc:478
void initialize() override
Fill in the options descriptions accordingly.
Definition: super.cc:49
COORD_TYPE coordtype_enum() const
Return the coordinate type recasted as the CASM defined enum.
Definition: Handlers.cc:398
const std::vector< std::string > & config_strs() const
Definition: Handlers.cc:382
const std::vector< std::string > & supercell_strs() const
Returns the list of the supercells for add_scelnames_suboption()
Definition: Handlers.cc:390
const std::string & unit_scel_str() const
Definition: super.cc:41
std::vector< fs::path > m_transf_mat_paths
Definition: Handlers.hh:786
const fs::path & struct_path() const
Definition: super.cc:39
double tolerance() const
Definition: super.cc:47
Index min_vol() const
Definition: super.cc:45
const std::vector< fs::path > & selection_paths() const
Returns the string corresponding to add_config_suboption()
Definition: Handlers.cc:342
const std::vector< fs::path > & transf_mat_paths() const
Definition: super.cc:35
ConfigMapperResult import_structure(SimpleStructure const &_struc, Configuration const *hint_ptr=nullptr, std::vector< DoFKey > const &_hint_dofs={ "occ"}) const
imports structure specified by '_struc' into primclex()
std::string name() const
Definition: Named_impl.hh:14
boost::iterator_range< iterator > selected()
Definition: Selection.cc:237
void error(const std::string &what)
Definition: Log.hh:129
PrimClex is the top-level data structure for a CASM project.
Definition: PrimClex.hh:55
Structure specifies the lattice and atomic basis of a crystal.
Definition: Structure.hh:30
const SymGroup & point_group() const
Definition: Structure.cc:110
const Lattice & lattice() const
Definition: Structure.hh:100
Represents a supercell of the primitive parent crystal structure.
Definition: Supercell.hh:51
const Lattice & lattice() const
The super lattice.
Definition: Supercell.cc:239
Eigen::Matrix3l transf_mat() const
Definition: Supercell.cc:235
std::pair< DB::DatabaseIterator< Supercell >, bool > insert() const
Insert the canonical form of this into the database.
Definition: Supercell.cc:283
SymGroup is a collection of symmetry operations that satisfy the group property The symmetry operatio...
Definition: SymGroup.hh:42
Print POSCAR with formating options.
Definition: VaspIO.hh:40
void set_coord_mode(COORD_TYPE mode)
Set coordinate mode.
Definition: VaspIO.hh:76
void set_atom_names_off()
Do not print atom names line.
Definition: VaspIO.hh:85
void set_atom_names_on()
Print atom names line.
Definition: VaspIO.hh:88
void sort()
Default sort is by atom name.
Definition: VaspIO.cc:50
void print(std::ostream &sout) const
Print POSCAR to stream.
Definition: VaspIO.cc:56
BasicStructure specifies the lattice and atomic basis of a crystal.
static BasicStructure from_poscar_stream(std::istream &poscar_stream, double tol=TOL)
const std::string & title() const
void set_title(std::string const &_title)
Set the title of the structure.
const Lattice & lattice() const
COORD_MODE specifies the current coordinate mode (Fractional or Cartesian)
const Eigen::Matrix3d & lat_column_mat() const
3x3 matrix with lattice vectors as its columne
Definition: Lattice.hh:110
double tol() const
Definition: Lattice.hh:195
PrimClex & make_primclex_if_not(const CommandArgs &args, std::unique_ptr< PrimClex > &uniq_primclex)
If !_primclex, construct new PrimClex stored in uniq_primclex, then return reference to existing or c...
Eigen::Matrix3i enforce_min_volume(const Lattice &unit, const Eigen::Matrix3i &T, const SymOpVector &point_grp, Index volume, bool fix_shape=false)
Return a transformation matrix that ensures a supercell of at least some volume.
Lattice make_superlattice(const Lattice &lat, const Eigen::Matrix< IntegralType, 3, 3, Options > &transf_mat)
Returns a super Lattice. Transformation matrix must be integer.
Definition: Lattice.hh:287
std::pair< bool, Eigen::Matrix3d > is_superlattice(const Lattice &scel, const Lattice &unit, double tol)
Check if scel is a superlattice of unitcell unit and some integer transformation matrix T.
Definition: Lattice.cc:836
std::pair< Eigen::MatrixXi, Eigen::MatrixXi > hermite_normal_form(const Eigen::MatrixXi &M)
Return the hermite normal form, M == H*V.
Eigen::CwiseUnaryOp< decltype(Local::iround_l< typename Derived::Scalar >), const Derived > iround(const Eigen::MatrixBase< Derived > &val)
Round Eigen::MatrixXd to Eigen::MatrixXi.
DB::Database< T > & db() const
Definition: PrimClex.cc:302
const PrimType & prim() const
const Access to primitive Structure
Definition: PrimClex.cc:262
void commit(ProjectSettings const &set)
std::string scelname(const Structure &prim, const Lattice &superlat)
Make supercell name name [deprecated].
Definition: Supercell.cc:497
ConfigIO::GenericConfigFormatter< jsonParser > config()
Definition: ConfigIO.cc:777
ConfigIO::GenericConfigFormatter< std::string > configname()
Constructs DataFormmaterDictionary containing all Configuration DatumFormatters.
Definition: ConfigIO.cc:563
Lattice equivalent(Lattice const &in_lat, SymOpVector const &point_grp, double compare_tol)
Lattice make_equivalent_superduperlattice(LatIterator begin, LatIterator end, SymOpIterator op_begin, SymOpIterator op_end)
returns Lattice that is smallest possible superlattice of all input Lattice
Definition: SymTools.hh:130
BasicStructure make_superstructure(const BasicStructure &tiling_unit, const Eigen::Matrix< IntegralType, 3, 3, Options > &transformation_matrix)
std::pair< OpIterator, Eigen::Matrix3d > is_equivalent_superlattice(const Object &scel, const Object &unit, OpIterator begin, OpIterator end, double tol)
Definition: SymTools.hh:110
Main CASM namespace.
Definition: APICommand.hh:8
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.
COORD_TYPE
Definition: enum.hh:6
Log & log()
Definition: Log.hh:424
Eigen::Matrix3d Matrix3d
const double TOL
Definition: definitions.hh:30
int super_command(const CommandArgs &args)
Definition: super.cc:129
Iterator find(Iterator begin, Iterator end, const T &value, BinaryCompare q)
Equivalent to std::find(begin, end, value), but with custom comparison.
Definition: algorithm.hh:16
INDEX_TYPE Index
For long integer indexing:
Definition: definitions.hh:39
Log & err_log()
Definition: Log.hh:426
PrimClex * primclex
Definition: settings.cc:135
Data structure holding basic CASM command info.
Struct with optional parameters for Config Mapping Specifies default parameters for all values,...