1 #ifndef CASM_Configuration
2 #define CASM_Configuration
26 using xtal::UnitCellCoord;
28 class IntegralCluster;
32 class ConfigIsEquivalent;
33 template <
typename ConfigType,
typename IsEqualImpl>
34 class GenericConfigCompare;
89 std::shared_ptr<Supercell const>
const &_supercell_ptr);
101 explicit Configuration(std::shared_ptr<Supercell const>
const &_supercell_ptr,
107 const std::shared_ptr<Supercell const> &_supercell_ptr);
112 const std::shared_ptr<Supercell const> &_supercell_ptr,
double _tol);
442 std::string
name()
const;
bool has_occupation() const
Eigen::VectorXi const & occupation() const
Const reference occupation values.
int & occ(Index i)
Reference occupation value on site i.
void set_occupation(Eigen::Ref< const Eigen::VectorXi > const &_occupation)
Set occupation values.
Class for comparison of Configurations (with the same Supercell)
Eigen::VectorXd num_each_component() const
bool is_canonical() const
const ConfigDoF & configdof() const
const Access the DoF
static Configuration zeros(const std::shared_ptr< Supercell const > &_supercell_ptr)
void set_occupation(Eigen::Ref< const Eigen::VectorXi > const &_occupation)
Set occupant variables.
Configuration primitive() const
Return the primitive Configuration.
Eigen::VectorXi const & occupation() const
Occupant variables.
std::shared_ptr< Supercell const > m_supercell_ptr
Used when constructing temporary Configuration in non-canonical Supercell.
std::vector< Eigen::VectorXd > sublattice_composition() const
Configuration & apply_sym(const PermuteIterator &it)
Transform Configuration from PermuteIterator via *this = permute_iterator * *this.
void clear()
Clear all DoF.
ConfigIsEquivalent equal_to() const
bool eq_impl(const Configuration &B) const
operator== comparison of Configuration, via ConfigEqual
bool has_occupation() const
True if Configuration has occupation DoF.
Configuration(std::shared_ptr< Supercell const > const &_supercell_ptr)
Eigen::VectorXi num_each_molecule() const
Eigen::VectorXd true_composition() const
const Lattice & ideal_lattice() const
std::vector< PermuteIterator > point_group() const
Returns the point group that leaves the Configuration unchanged.
PermuteIterator find_translation() const
Returns a PermuteIterator corresponding to the first non-zero pure translation that maps the Configur...
std::vector< Eigen::VectorXi > sublat_num_each_molecule() const
std::vector< PermuteIterator > invariant_subgroup() const
Returns the subgroup of the Supercell factor group that leaves the Configuration unchanged.
Index linear_index(const UnitCellCoord &bijk) const
Return the linear index corresponding to integral coordinates.
ConfigInsertResult insert(bool primitive_only=false) const
Insert this configuration (in primitive & canonical form) in the database.
Supercell const * m_supercell
Pointer to the Supercell for this Configuration.
std::string generate_name_impl() const
Returns a Configuration name.
UnitCellCoord uccoord(Index site_l) const
Get the UnitCellCoord for a given linear site index.
const Molecule & mol(Index site_l) const
Molecule on site l.
Eigen::VectorXd param_composition() const
Returns parametric composition, as calculated using PrimClex::param_comp.
const int & occ(Index site_l) const
Occupant variable on site l.
void init_occupation()
Set occupant variables to background structure.
static std::pair< std::string, std::string > split_name(std::string configname)
Split configuration name string into scelname and config id.
int multiplicity() const
Get symmetric multiplicity, excluding translations.
ConfigDoF & configdof()
Access the DoF.
ConfigCompare less() const
const Supercell & supercell() const
Get the Supercell for this Configuration.
ConfigDoF m_configdof
Degrees of Freedom.
Eigen::VectorXd composition() const
Configuration in_canonical_supercell() const
Returns the canonical form Configuration in the canonical Supercell.
void set_occ(Index site_l, int val)
Set occupant variable on site l.
std::string point_group_name() const
Returns the point group that leaves the Configuration unchanged.
int sublat(Index site_l) const
Get the basis site index for a given linear linear site index.
Index size() const
Returns number of sites, NOT the number of primitives that fit in here.
bool is_primitive() const
Check if this is a primitive Configuration.
bool operator<(const Configuration &B) const
Clear occupation.
std::vector< PermuteIterator > factor_group() const
Returns the subgroup of the Supercell factor group that leaves the Configuration unchanged.
Class for less than comparison of Configurations implemented via a ConfigTypeIsEqual class that also ...
Implements PrimClex dependent functions.
PrimClex is the top-level data structure for a CASM project.
Represents a supercell of the primitive parent crystal structure.
SymOp is the Coordinate representation of a symmetry operation it keeps fraction (FRAC) and Cartesian...
Class representing a Molecule.
Base class for CRTP pattern.
double reference_energy(const Configuration &config)
Returns the reference energy, normalized per unit cell.
bool has_relaxed_magmom(const Configuration &_config)
Eigen::MatrixXd relaxed_forces(const Configuration &_config)
relaxed forces of configuration, determined from DFT (eV/Angstr.), as a 3xN matrix
bool has_volume_relaxation(const Configuration &_config)
double formation_energy(const Configuration &config)
Returns the formation energy, normalized per unit cell.
bool has_reference_energy(const Configuration &_config)
double formation_energy_per_species(const Configuration &config)
Returns the formation energy, normalized per species.
double clex_formation_energy_per_species(const Configuration &config)
Returns the formation energy, normalized per species.
double atomic_deformation(const Configuration &_config)
Cost function that describes the degree to which basis sites have relaxed.
Eigen::VectorXi num_each_molecule(const ConfigDoF &configdof, const Supercell &scel)
Returns num_each_molecule(molecule_type), where 'molecule_type' is ordered as Structure::get_struc_mo...
ConfigCanonicalForm< HasSupercell< Comparisons< Calculable< CRTPBase< Configuration > > > > > ConfigurationBase
std::string diff_trans_endpoint_of(const Configuration &_config)
returns which diff_trans _config is an endpoint of
Eigen::VectorXd correlations(const Configuration &config, Clexulator const &clexulator)
Returns correlations using 'clexulator'.
double clex_formation_energy(const Configuration &config)
Returns the formation energy, normalized per unit cell.
Eigen::VectorXd corr_contribution(Index linear_unitcell_index, const Configuration &config, Clexulator const &clexulator)
Returns correlation contribution from a single unit cell, not normalized.
bool is_diff_trans_endpoint(const Configuration &_config)
returns true if _config is an endpoint of any existing diff_trans_config
Eigen::VectorXd comp(const Configuration &config)
Returns parametric composition, as calculated using PrimClex::param_comp.
Configuration closest_setting(const Configuration &_config1, const Configuration &_config2)
Configuration sub_configuration(std::shared_ptr< Supercell const > sub_scel_ptr, const Configuration &super_config, const UnitCell &origin=UnitCell(0, 0, 0))
Returns the sub-configuration that fills a particular Supercell.
double rms_force(const Configuration &_config)
Root-mean-square forces of relaxed configurations, determined from DFT (eV/Angstr....
bool has_lattice_deformation(const Configuration &_config)
bool has_rms_force(const Configuration &_config)
double n_species(const Configuration &config)
Returns the total number species per unit cell.
bool is_primitive(const Configuration &_config)
returns true if _config describes primitive cell of the configuration it describes
double volume_relaxation(const Configuration &_config)
Change in volume due to relaxation, expressed as the ratio V/V_0.
double reference_energy_per_species(const Configuration &config)
Returns the reference energy, normalized per species.
Eigen::MatrixXd gradcorrelations(const Configuration &config, Clexulator const &clexulator, DoFKey &key)
Returns gradient correlations using 'clexulator', with respect to DoF 'dof_type'.
Eigen::VectorXd comp_n(const Configuration &config)
Returns the composition, as number of each species per unit cell.
IntegralCluster config_diff(const Configuration &_config1, const Configuration &_config2)
Returns an Integral Cluster representing the perturbed sites between the configs.
double relaxed_magmom(const Configuration &_config)
Returns the relaxed magnetic moment, normalized per unit cell.
Eigen::VectorXd site_frac(const Configuration &config)
Returns the composition as site fraction, in the order of Structure::get_struc_molecule.
Configuration make_configuration(const PrimClex &primclex, std::string name)
Make Configuration from name string.
double lattice_deformation(const Configuration &_config)
Cost function that describes the degree to which lattice has relaxed.
double energy(const Configuration &config)
Returns the energy, normalized per unit cell.
double n_vacancy(const Configuration &config)
Returns the vacancy composition, as number per unit cell.
Configuration config_clip(const Configuration &_config, const Configuration &_bg, IntegralCluster &_clust)
Returns a Configuration with the sites in _clust clipped from _config and placed in _bg.
Eigen::VectorXd point_corr(Index linear_unitcell_index, Index neighbor_index, const Configuration &config, Clexulator const &clexulator)
Returns point correlations from a single site, normalized by cluster orbit size.
double energy_per_species(const Configuration &config)
Returns the energy, normalized per species.
bool has_formation_energy(const Configuration &_config)
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...
bool has_atomic_deformation(const Configuration &_config)
bool is_canonical(const Configuration &_config)
returns true if no symmetry transformation applied to _config will increase its lexicographic order
bool has_relaxed_mag_basis(const Configuration &_config)
double relaxed_magmom_per_species(const Configuration &_config)
Returns the relaxed magnetic moment, normalized per species.
bool has_energy(const Configuration &_config)
ConfigIO::GenericConfigFormatter< jsonParser > config()
ConfigIO::GenericConfigFormatter< std::string > configname()
Constructs DataFormmaterDictionary containing all Configuration DatumFormatters.
GenericDatumFormatter< std::string, DataObject > name()
INDEX_TYPE Index
For long integer indexing:
Holds results of Configuration::insert.
bool insert_canonical
True if canonical configuration did not exist before insertion.
iterator primitive_it
Iterator pointing at primitive.
bool insert_primitive
True if primitive did not exist before insertion.
iterator canonical_it
Iterator pointing at canonical, if existing.
DB::DatabaseIterator< Configuration > iterator
Operations that transform a canonical primitive configuration to any equivalent.
Configuration prim_canon_config
PermuteIterator from_canonical_config
Eigen::Matrix3i transf_mat
RefToCanonicalPrim(const Configuration &_config)
Get operations that transform canonical primitive to this.
Implements other comparisons in terms of '<'.