CASM
AClustersApproachtoStatisticalMechanics
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules
Configuration.hh
Go to the documentation of this file.
1 #ifndef CASM_Configuration
2 #define CASM_Configuration
3 
4 #include <map>
5 
6 #include "casm/external/boost.hh"
7 
10 #include "casm/container/Array.hh"
13 #include "casm/clex/Properties.hh"
14 #include "casm/clex/Correlation.hh"
16 #include "casm/clex/ConfigDoF.hh"
18 
19 namespace CASM {
20 
21  class PrimClex;
22  class Supercell;
23  class UnitCellCoord;
24  class Clexulator;
25  class FillSupercell;
26  template<typename ConfigType, typename PrimClexType>
27  class ConfigIterator;
28 
29 
30  namespace ConfigIO {
31  class Selected;
32  }
33  template <bool IsConst> class ConfigSelection;
34  typedef ConfigSelection<true> ConstConfigSelection;
35 
36  template<typename DataObject> struct QueryTraits;
37 
38  template<>
40  static const std::string name;
42  typedef ConstConfigSelection Selection;
43  };
44 
52 
53 
62 
64 
67 
69  config_const_iterator primitive_it;
70 
73 
75  config_const_iterator canonical_it;
76 
77  };
78 
81  class Configuration : public Comparisons<Configuration> {
82  private:
83 
86 
89 
90 
92 
93  // Configuration id is the index into Supercell::config_list
94  std::string id;
95 
98 
100  std::shared_ptr<Supercell> m_supercell_ptr;
101 
105 
106 
107  // symmetric multiplicity (i.e., size of configuration's factor_group)
109 
110 
112 
113  // 'occupation' is a list of the indices describing the occupants in each crystal site.
114  // get_prim().basis[ get_b(i) ].site_occupant[ occupation[i]] -> Molecule on site i
115  // This means that for the background structure, 'occupation' is all 0
116 
117  // Configuration sites are arranged by basis, and then prim:
118  // occupation: [basis0 |basis1 |basis2 |...] up to prim.basis.size()
119  // basis0: [prim0|prim1|prim2|...] up to supercell.volume()
120  //
122 
123 
127  /* calculated:
128  * calculated["energy"]
129  * calculated["relaxed_energy"]
130  *
131  * generated:
132  * generated["is_groundstate"]
133  * generated["dist_from_hull"]
134  * generated["sublat_struct_fact"]
135  * generated["struct_fact"]
136  */
138  Properties calculated; //Stuff you got directly from your DFT calculations
139  Properties generated; //Everything else you came up with through casm
140 
141 
143 
146 
148  mutable std::string m_name;
149 
150  public:
154 
155  //********* CONSTRUCTORS *********
156 
158  Configuration(Supercell &_supercell, const jsonParser &source = jsonParser(), const ConfigDoF &_dof = ConfigDoF());
159 
161  Configuration(const jsonParser &json, Supercell &_supercell, Index _id);
162 
163 
165  //Configuration(Supercell &_supercell, std::string con_name, bool select, const jsonParser &source = jsonParser());
166 
167 
168  //********** DESTRUCTORS *********
169 
170  //********** MUTATORS ***********
171 
173  void set_multiplicity(int m) {
174  multiplicity = m;
175  }
176 
177  void set_id(Index _id);
178 
179  void set_source(const jsonParser &source);
180 
181  void push_back_source(const jsonParser &source);
182 
183  // ** Degrees of Freedom **
184  //
185  // ** Note: Properties and correlations are not automatically updated when dof are changed, **
186  // ** nor are the written records automatically updated **
187 
191  void clear();
192 
193 
197  void init_occupation();
198 
214  void set_occupation(const Array<int> &newoccupation);
215 
228  void set_occ(Index site_l, int val);
229 
233  void clear_occupation();
234 
235 
239  void init_displacement();
240 
251  void set_displacement(const displacement_matrix_t &_disp);
252 
259  void set_disp(Index site_l, const Eigen::VectorXd &_disp);
260 
264  void clear_displacement();
265 
266 
270  void init_deformation();
271 
284  void set_deformation(const Eigen::Matrix3d &_deformation);
285 
289  void clear_deformation();
290 
291 
293  bool is_primitive() const;
294 
298 
300  Configuration primitive() const;
301 
302 
304  bool is_canonical() const;
305 
308 
311 
314 
317 
319  ConfigInsertResult insert(bool primitive_only) const;
320 
323  std::vector<PermuteIterator> factor_group() const;
324 
326  SymGroup point_group() const;
327 
329  Configuration fill_supercell(Supercell &scel, const SymOp &op) const;
330 
332  Configuration fill_supercell(Supercell &scel, const SymGroup &g) const;
333 
334 
335  // ** Properties **
336  //
337  // ** Note: DeltaProperties are automatically updated, but not written upon changes **
338 
340  //void read_calculated();
341 
342  void set_calc_properties(const jsonParser &json);
343 
344  bool read_calc_properties(jsonParser &parsed_props) const;
345 
346  void set_selected(bool _selected) {
347  m_selected = _selected;
348  }
349 
350 
351  //********** ACCESSORS ***********
352 
353  const Lattice &ideal_lattice()const;
354 
355  std::string get_id() const;
356 
360  int get_multiplicity()const {
361  return multiplicity;
362  }
363 
365  std::string name() const;
366 
367  std::string calc_status() const;
368 
369  std::string failure_type() const;
370 
371  const jsonParser &source() const;
372 
373  fs::path get_path() const;
374 
376  Index size() const;
377 
379  const Structure &get_prim() const;
380 
382  bool selected() const {
383  return m_selected;
384  }
385 
387  PrimClex &get_primclex() const;
388 
390  Supercell &get_supercell() const;
391 
393  double crystallography_tol() const;
394 
396  UnitCellCoord get_uccoord(Index site_l) const;
397 
399  int get_b(Index site_l) const;
400 
402  const ConfigDoF &configdof() const {
403  return m_configdof;
404  }
405 
410  _invalidate_id();
411  return m_configdof;
412  }
413 
415  bool has_occupation() const {
416  return configdof().has_occupation();
417  }
418 
431  const Array<int> &occupation() const {
432  return configdof().occupation();
433  }
434 
445  const int &occ(Index site_l) const {
446  return configdof().occ(site_l);
447  }
448 
456  const Molecule &get_mol(Index site_l) const;
457 
458 
460  bool has_displacement() const {
461  return configdof().has_displacement();
462  }
463 
471  const displacement_matrix_t &displacement() const {
472  return configdof().displacement();
473  }
474 
480  const_displacement_t disp(Index site_l) const {
481  return configdof().disp(site_l);
482  }
483 
496  const Eigen::Matrix3d &deformation() const {
497  return configdof().deformation();
498  }
499 
501  bool has_deformation() const {
502  return configdof().has_deformation();
503  }
504 
505 
506  //fs::path get_reference_state_dir() const;
507 
508  //const Properties &ref_properties() const;
509 
510  const Properties &calc_properties() const;
511 
512  //const DeltaProperties &delta_properties() const;
513 
514  const Properties &generated_properties() const;
515 
516 
517  //const Correlation &get_correlations() const;
518 
519 
520  // Returns composition on each sublattice: sublat_comp[ prim basis site / sublattice][ molecule_type]
521  // molucule_type is ordered as in the Prim structure's site_occupant list for that basis site (includes vacancies)
523 
524  // Returns number of each molecule by sublattice:
525  // sublat_num_each_molecule[ prim basis site / sublattice ][ molecule_type]
526  // molucule_type is ordered as in the Prim structure's site_occupant list for that basis site
528 
529  // Returns composition, not counting vacancies
530  // composition[ molecule_type ]: molecule_type ordered as prim structure's get_struc_molecule(), with [Va]=0.0
532 
533  // Returns composition, including vacancies
534  // composition[ molecule_type ]: molecule_type ordered as prim structure's get_struc_molecule()
536 
539 
542 
546 
547  //-----------------------------------
548  //Structure Factor
550  Eigen::VectorXd get_struct_fact_intensities(const Eigen::VectorXd &component_intensities) const;
551 
553  void calc_struct_fact();
554  void calc_sublat_struct_fact(const Eigen::VectorXd &intensities);
555  void calc_struct_fact(const Eigen::VectorXd &intensities);
556 
557  Eigen::MatrixXcd sublat_struct_fact();
559 
560  //********* IO ************
561 
565  jsonParser &write(jsonParser &json) const;
566 
568  void write_pos() const;
569 
570  // Va_mode description
571  // 0 print no information about the vacancies
572  // 1 print only the coordinates of the vacancies
573  // 2 print the number of vacancies and the coordinates of the vacancies
574  void print(std::ostream &stream, COORD_TYPE mode, int Va_mode = 0, char term = '\n', int prec = 10, int pad = 5) const;
575 
576  void print_occupation(std::ostream &stream) const;
577 
578  void print_config_list(std::ostream &stream, int composition_flag) const;
579 
580  void print_composition(std::ostream &stream) const;
581 
582  void print_true_composition(std::ostream &stream) const;
583 
584  void print_sublattice_composition(std::ostream &stream) const;
585 
586 
587  fs::path calc_dir() const;
589  fs::path calc_status_path() const;
591  fs::path get_pos_path() const;
592 
601  bool is_equivalent(const Configuration &B) const;
602 
606  bool operator<(const Configuration &B) const;
607 
609  static std::pair<std::string, Index> split_name(std::string configname);
610 
611  private:
613  int &_occ(Index site_l) {
614  return m_configdof.occ(site_l);
615  }
616 
617  displacement_t _disp(Index site_l) {
618  return m_configdof.disp(site_l);
619  }
620 
621  void _invalidate_id() {
622  id = "none";
623  m_name.clear();
624  }
625 
626  void _generate_name() const;
627 
629 
634  bool _eq(const Configuration &B) const;
635 
644  void read(const jsonParser &json);
645 
647  void read_dof(const jsonParser &json);
648  void read_properties(const jsonParser &json);
649 
651  jsonParser &write_dof(jsonParser &json) const;
652  jsonParser &write_source(jsonParser &json) const;
653  jsonParser &write_pos(jsonParser &json) const;
655  jsonParser &write_properties(jsonParser &json) const;
656 
657  //bool reference_states_exist() const;
658  //void read_reference_states(Array<Properties> &ref_state_prop, Array<Eigen::VectorXd> &ref_state_comp) const;
659  //void generate_reference_scalar(std::string propname, const Array<Properties> &ref_state_prop, const Array<Eigen::VectorXd> &ref_state_comp);
660 
661  };
662 
663  // Calculate transformed ConfigDoF from PermuteIterator via
664  // apply(permute_iterator, dof)
665  Configuration &apply(const PermuteIterator &it, Configuration &config);
666 
668  const Configuration &super_config,
669  const UnitCell &origin = UnitCell(0, 0, 0));
670 
673 
675  Correlation correlations(const Configuration &config, Clexulator &clexulator);
676 
678  Eigen::VectorXd comp(const Configuration &config);
679 
681  Eigen::VectorXd comp_n(const Configuration &config);
682 
684  double n_vacancy(const Configuration &config);
685 
687  double n_species(const Configuration &config);
688 
692 
694  Eigen::VectorXd site_frac(const Configuration &config);
695 
697  double relaxed_energy(const Configuration &config);
698 
700  double relaxed_energy_per_species(const Configuration &config);
701 
703  double reference_energy(const Configuration &config);
704 
706  double reference_energy_per_species(const Configuration &config);
707 
709  double formation_energy(const Configuration &config);
710 
712  double formation_energy_per_species(const Configuration &config);
713 
715  double clex_formation_energy(const Configuration &config);
716 
718  double clex_formation_energy_per_species(const Configuration &config);
719 
721  bool is_calculated(const Configuration &config);
722 
724  double rms_force(const Configuration &_config);
725 
727  double basis_deformation(const Configuration &_config);
728 
730  double lattice_deformation(const Configuration &_config);
731 
733  double volume_relaxation(const Configuration &_config);
734 
736  double relaxed_magmom(const Configuration &_config);
737 
739  double relaxed_magmom_per_species(const Configuration &_config);
740 
743  /* std::vector<double> relaxed_mag_basis(const Configuration &_config); */
744 
746  Eigen::VectorXd relaxed_mag(const Configuration &_config);
747 
748 
749 
751  bool is_primitive(const Configuration &_config);
752 
754  bool is_canonical(const Configuration &_config);
755 
757  inline
758  std::string calc_status(const Configuration &_config) {
759  return _config.calc_status();
760  }
761 
762  // \brief Reason for calculation failure.
763  inline
764  std::string failure_type(const Configuration &_config) {
765  return _config.failure_type();
766  }
767 
768  bool has_relaxed_energy(const Configuration &_config);
769 
770  bool has_reference_energy(const Configuration &_config);
771 
772  bool has_formation_energy(const Configuration &_config);
773 
774  bool has_rms_force(const Configuration &_config);
775 
776  bool has_basis_deformation(const Configuration &_config);
777 
778  bool has_lattice_deformation(const Configuration &_config);
779 
780  bool has_volume_relaxation(const Configuration &_config);
781 
782  bool has_relaxed_magmom(const Configuration &_config);
783 
784  bool has_relaxed_mag_basis(const Configuration &_config);
785 
786  inline
787  bool has_calc_status(const Configuration &_config) {
788  return !_config.calc_status().empty();
789  }
790 
791  inline
792  bool has_failure_type(const Configuration &_config) {
793  return !_config.failure_type().empty();
794  }
795 
797 
798  public:
799 
801  FillSupercell(Supercell &_scel, const SymOp &_op);
802 
805  FillSupercell(Supercell &_scel, const Configuration &motif, double tol);
806 
807  Configuration operator()(const Configuration &motif) const;
808 
811  const SymOp *find_symop(const Configuration &motif, double tol);
812 
813  const SymOp &symop() const {
814  return *m_op;
815  }
816 
817 
818  private:
819 
820  void _init(Supercell &_motif_scel) const;
821 
823  const SymOp *m_op;
824 
826  mutable std::vector<std::vector<Index> > m_index_table;
827 
828  };
829 
830  std::ostream &operator<<(std::ostream &sout, const Configuration &c);
831 
832 
833  inline
835  _config.set_calc_properties(jsonParser());
836  }
837 
840 }
841 
842 #endif
Eigen::MatrixXd MatrixXd
void clear_deformation()
Clear applied strain.
void clear_displacement()
Clear displacement.
bool insert_primitive
True if primitive did not exist before insertion.
Eigen::MatrixXd pad(const Eigen::MatrixXd &M, int n)
Construct a matrix consisting of blocks M and Identity(n,n)
Definition: PCA.hh:88
Configuration fill_supercell(Supercell &scel, const SymOp &op) const
Fills supercell 'scel' with reoriented configuration, as if by apply(op,*this)
ConfigInsertResult insert(bool primitive_only) const
Insert this in the canonical Supercell.
double rms_force(const Configuration &_config)
Root-mean-square forces of relaxed configurations, determined from DFT (eV/Angstr.)
Eigen::VectorXd get_param_composition() const
Returns parametric composition, as calculated using PrimClex::param_comp.
const_displacement_t disp(Index site_l) const
Occupant displacement.
bool selected() const
True if this Configuration is currently selected in the MASTER config list.
static std::pair< std::string, Index > split_name(std::string configname)
Split configuration name string into scelname and config index.
int & occ(Index i)
Definition: ConfigDoF.hh:61
void print(std::ostream &stream, COORD_TYPE mode, int Va_mode=0, char term= '\n', int prec=10, int pad=5) const
void reset_properties(ConfigDoF &_dof)
Definition: ConfigDoF.hh:202
Eigen::MatrixXcd sublat_struct_fact()
void init_occupation()
Set occupant variables to background structure.
static const std::string name
UnitCellCoord get_uccoord(Index site_l) const
Get the UnitCellCoord for a given linear site index.
Eigen::VectorXd relaxed_mag(const Configuration &_config)
Returns the relaxed magnetic moment for each molecule.
const int & occ(Index site_l) const
Occupant variable on site l.
Supercell * supercell
const pointer to the (non-const) Supercell for this Configuration
Implements other comparisons in terms of '<'.
Definition: Comparisons.hh:23
void set_deformation(const Eigen::Matrix3d &_deformation)
Set applied strain.
PrimClex & get_primclex() const
Get the PrimClex for this Configuration.
double relaxed_magmom(const Configuration &_config)
Returns the relaxed magnetic moment, normalized per unit cell.
void read_dof(const jsonParser &json)
Functions used to perform read()
void read_properties(const jsonParser &json)
const Molecule & get_mol(Index site_l) const
Molecule on site l.
std::string calc_status() const
bool has_basis_deformation(const Configuration &_config)
ConfigIterator< const Configuration, const PrimClex > config_const_iterator
bool is_equivalent(const Configuration &B) const
Check if Configuration are equivalent wrt the prim's factor group.
const Array< int > & occupation() const
Occupant variables.
bool has_formation_energy(const Configuration &_config)
Structure specifies the lattice and atomic basis of a crystal.
Definition: Structure.hh:29
Configuration make_configuration(PrimClex &primclex, std::string name)
Make Configuration from name string.
PrimClex * primclex
Definition: settings.cc:101
fs::path calc_properties_path() const
config_const_iterator primitive_it
Iterator pointing at primitive.
const SymOp & symop() const
Configuration canonical_form() const
Returns the canonical form Configuration in the same Supercell.
void print_sublattice_composition(std::ostream &stream) const
Eigen::MatrixXd struct_fact()
void print_occupation(std::ostream &stream) const
Unit Cell Coordinates.
bool has_lattice_deformation(const Configuration &_config)
Correlation correlations(const ConfigDoF &configdof, const Supercell &scel, Clexulator &clexulator)
Returns correlations using 'clexulator'. Supercell needs a correctly populated neighbor list...
Definition: ConfigDoF.cc:200
double clex_formation_energy_per_species(const Configuration &config)
Returns the formation energy, normalized per species.
bool has_relaxed_energy(const Configuration &_config)
const jsonParser & source() const
const ConfigDoF & configdof() const
const Access the DoF
void push_back_source(const jsonParser &source)
config_const_iterator canonical_it
Iterator pointing at canonical, if existing.
Eigen::VectorXd get_num_each_component() const
ConfigDoF m_configdof
Degrees of Freedom.
void print_true_composition(std::ostream &stream) const
void set_displacement(const displacement_matrix_t &_disp)
Set occupant displacements.
double n_species(const Configuration &config)
Returns the total number species per unit cell.
int get_b(Index site_l) const
Get the basis site index for a given linear linear site index.
fs::path calc_status_path() const
Configuration primitive() const
Return the primitive Configuration.
Main CASM namespace.
Definition: complete.cpp:8
void set_disp(Index site_l, const Eigen::VectorXd &_disp)
Set occupant displacements.
Eigen::VectorXd comp(const Configuration &config)
Returns parametric composition, as calculated using PrimClex::param_comp.
const Lattice & ideal_lattice() const
Represents a supercell of the primitive parent crystal structure.
Definition: Supercell.hh:37
jsonParser m_source
a jsonParser object indicating where this Configuration came from
Unit Cell Indices.
std::string m_name
Remember name.
const Properties & generated_properties() const
ReturnArray< Array< int > > get_sublat_num_each_molecule() const
displacement_t _disp(Index site_l)
fs::path calc_dir() const
const Eigen::Matrix3d & deformation() const
Definition: ConfigDoF.hh:85
Configuration operator()(const Configuration &motif) const
bool read_calc_properties(jsonParser &parsed_props) const
double relaxed_energy_per_species(const Configuration &config)
Returns the relaxed energy, normalized per species.
ConfigDoF & configdof()
Access the DoF.
Eigen::VectorXd get_struct_fact_intensities() const
bool insert_canonical
True if canonical configuration did not exist before insertion.
const Structure & get_prim() const
Get the primitive Structure for this Configuration.
SymGroup is a collection of symmetry operations that satisfy the group property The symmetry operatio...
Definition: SymGroup.hh:33
std::string get_id() const
bool _eq(const Configuration &B) const
Equality comparison of Configuration, via ConfigEqual.
Supercell * m_motif_scel
Eigen::VectorXd comp_n(const ConfigDoF &configdof, const Supercell &scel)
Returns comp_n, the number of each molecule per primitive cell, ordered as Structure::get_struc_molec...
Definition: ConfigDoF.cc:312
std::shared_ptr< Supercell > m_supercell_ptr
double tol
void _init(Supercell &_motif_scel) const
bool has_occupation() const
True if Configuration has occupation DoF.
displacement_matrix_t::ColXpr displacement_t
Definition: ConfigDoF.hh:35
double formation_energy(const Configuration &config)
Returns the formation energy, normalized per unit cell.
double clex_formation_energy(const Configuration &config)
Returns the formation energy, normalized per unit cell.
bool is_canonical() const
Check if Configuration is in the canonical form.
const Properties & calc_properties() const
std::ostream & operator<<(std::ostream &_stream, const FormattedPrintable &_formatted)
SymOp is the Coordinate representation of a symmetry operation it keeps fraction (FRAC) and Cartesian...
Definition: SymOp.hh:28
void print_composition(std::ostream &stream) const
displacement_matrix_t::ConstColXpr const_displacement_t
Definition: ConfigDoF.hh:36
fs::path get_pos_path() const
Path to various files.
std::string failure_type(const Configuration &_config)
std::vector< std::vector< Index > > m_index_table
bool has_relaxed_mag_basis(const Configuration &_config)
EigenIndex Index
For long integer indexing:
const SymOp * m_op
Evaluates correlations.
Definition: Clexulator.hh:240
bool has_rms_force(const Configuration &_config)
void set_occupation(const Array< int > &newoccupation)
Set occupant variables.
bool has_deformation() const
True if Configuration has strain DoF.
std::string id
Identification.
A container class for the different degrees of freedom a Configuration might have.
Definition: ConfigDoF.hh:27
ReturnArray< int > get_num_each_molecule() const
Returns num_each_molecule[ molecule_type], where 'molecule_type' is ordered as Structure::get_struc_m...
bool is_canonical(const Configuration &_config)
returns true if _config no symmetry transformation applied to _config will increase its lexicographic...
double volume_relaxation(const Configuration &_config)
Change in volume due to relaxation, expressed as the ratio V/V_0.
bool is_calculated(const Configuration &config)
Return true if all current properties have been been calculated for the configuration.
fs::path get_path() const
const SymOp * find_symop(const Configuration &motif, double tol)
Find first SymOp in the prim factor group such that apply(op, motif) can be used to fill the Supercel...
Eigen::VectorXd VectorXd
void clear_occupation()
Clear occupation.
PrimClex is the top-level data structure for a CASM project.
Definition: PrimClex.hh:52
Holds results of Configuration::insert.
void set_id(Index _id)
void _generate_name() const
PermuteIterator find_translation() const
Returns a PermuteIterator corresponding to the first non-zero pure translation that maps the Configur...
std::string name() const
SCELV_A_B_C_D_E_F/i.
double relaxed_magmom_per_species(const Configuration &_config)
Returns the relaxed magnetic moment, normalized per species.
jsonParser & write(jsonParser &json) const
Index size() const
Returns number of sites, NOT the number of primitives that fit in here.
ReturnArray< Array< double > > get_sublattice_composition() const
jsonParser & write_param_composition(jsonParser &json) const
notstd::cloneable_ptr< FillSupercell > m_fill_canonical
Remember how to copy into the canonical Supercell.
Eigen::VectorXd relaxed_mag_basis(const Configuration &_config)
Returns the relaxed magnetic moment of each basis site.
Configuration(Supercell &_supercell, const jsonParser &source=jsonParser(), const ConfigDoF &_dof=ConfigDoF())
Construct a default Configuration.
ConfigIO::GenericConfigFormatter< std::string > configname()
Constructs DataFormmaterDictionary containing all Configuration DatumFormatters.
Definition: ConfigIO.cc:340
FillSupercell(Supercell &_scel, const SymOp &_op)
Constructor.
Eigen::VectorXd Correlation
Definition: Correlation.hh:9
void set_calc_properties(const jsonParser &json)
Read calculation results into the configuration.
jsonParser & write_properties(jsonParser &json) const
std::vector< PermuteIterator > factor_group() const
Returns the subgroup of the Supercell factor group that leaves the Configuration unchanged.
Eigen::VectorXd site_frac(const Configuration &config)
Returns the composition as site fraction, in the order of Structure::get_struc_molecule.
void write_pos() const
Write the POS file to get_pos_path.
double formation_energy_per_species(const Configuration &config)
Returns the formation energy, normalized per species.
double lattice_deformation(const Configuration &_config)
Cost function that describes the degree to which lattice has relaxed.
ConfigDoF::displacement_t displacement_t
Eigen::Matrix3d Matrix3d
ReturnArray< double > get_true_composition() const
PermuteIterator to_canonical() const
Returns the operation that applied to *this returns the canonical form.
Configuration in_canonical_supercell() const
Returns the canonical form Configuration in the canonical Supercell.
bool has_failure_type(const Configuration &_config)
void set_source(const jsonParser &source)
const displacement_matrix_t & displacement() const
Occupant displacements.
const displacement_matrix_t & displacement() const
Definition: ConfigDoF.hh:81
const Array< int > & occupation() const
Definition: ConfigDoF.hh:69
bool has_volume_relaxation(const Configuration &_config)
SymGroup point_group() const
Returns the point group that leaves the Configuration unchanged.
void print_config_list(std::ostream &stream, int composition_flag) const
void set_occ(Index site_l, int val)
Set occupant variable on site l.
double relaxed_energy(const Configuration &config)
Returns the relaxed energy, normalized per unit cell.
const Eigen::Matrix3d & deformation() const
Applied strain.
double basis_deformation(const Configuration &_config)
Cost function that describes the degree to which basis sites have relaxed.
void clear()
Clear all DoF.
ConfigDoF::displacement_matrix_t displacement_matrix_t
double reference_energy(const Configuration &config)
Returns the reference energy, normalized per unit cell.
bool has_deformation() const
Definition: ConfigDoF.hh:98
bool is_primitive() const
Check if this is a primitive Configuration.
ConfigDoF::const_displacement_t const_displacement_t
ConfigSelection< true > ConstConfigSelection
bool has_relaxed_magmom(const Configuration &_config)
void set_selected(bool _selected)
void set_multiplicity(int m)
Construct a Configuration with occupation specified by string 'con_name'.
A 'cloneable_ptr' can be used in place of 'unique_ptr'.
int get_multiplicity() const
Get symmetric multiplicity (i.e., size of configuration's factor_group)
double reference_energy_per_species(const Configuration &config)
Returns the reference energy, normalized per species.
Eigen::MatrixXd displacement_matrix_t
Definition: ConfigDoF.hh:32
int & _occ(Index site_l)
Convenience accessors:
Configuration sub_configuration(Supercell &sub_scel, const Configuration &super_config, const UnitCell &origin=UnitCell(0, 0, 0))
Returns the sub-configuration that fills a particular Supercell.
Supercell & get_supercell() const
Get the Supercell for this Configuration.
displacement_t disp(Index i)
Definition: ConfigDoF.hh:73
bool has_displacement() const
Definition: ConfigDoF.hh:102
std::string calc_status(const Configuration &_config)
Status of calculation.
ReturnArray< double > get_composition() const
jsonParser & write_dof(jsonParser &json) const
Functions used to perform write to config_list.json:
void read(const jsonParser &json)
Private members:
void init_displacement()
Set all occupant displacements to (0.,0.,0.)
double n_vacancy(const Configuration &config)
Returns the vacancy composition, as number per unit cell.
bool has_displacement() const
True if Configuration has displacement DoF.
PermuteIterator from_canonical() const
Returns the operation that applied to the the canonical form returns *this.
bool is_primitive(const Configuration &_config)
returns true if _config describes primitive cell of the configuration it describes ...
Object & apply(const Transform &f, Object &obj, Args &&...args)
bool operator<(const Configuration &B) const
Compare Configuration, via ConfigCompare.
std::string failure_type() const
Eigen::VectorXd species_frac(const Configuration &config)
Returns the composition as species fraction, with [Va] = 0.0, in the order of Structure::get_struc_mo...
jsonParser & write_source(jsonParser &json) const
A Configuration represents the values of all degrees of freedom in a Supercell.
bool has_reference_energy(const Configuration &_config)
Returns true if configuration is specified in given selection (default: MASTER)
void init_deformation()
Set applied strain to Eigen::Matrix3d::Zero()
bool has_occupation() const
Definition: ConfigDoF.hh:106
bool has_calc_status(const Configuration &_config)
double crystallography_tol() const
Get the PrimClex crystallography_tol.