CASM  1.1.0
A Clusters Approach to Statistical Mechanics
CASM::xtal::StrucMapper Class Reference

#include <StrucMapping.hh>

Detailed Description

A class for mapping an arbitrary 'child' crystal structure as a deformation of a 'parent' crystal structure StrucMapper manages options for the mapping algorithm and mapping cost function It also caches some information about supercell lattices to improve speed of mapping multiple children crystals onto a single parent structure

Definition at line 404 of file StrucMapping.hh.

Public Types

enum  Options {
  none = 0 , strict = (1u << 0) , robust = (1u << 1) , sym_strain = (1u << 2) ,
  sym_basis = (1u << 3) , soft_va_limit = (1u << 4)
}
 
using LatMapType = std::map< Index, std::vector< Lattice > >
 

Public Member Functions

 StrucMapper (StrucMapCalculatorInterface const &_calculator, double _lattice_weight=0.5, double _max_volume_change=0.5, int _options=0, double _cost_tol=TOL, double _min_va_frac=0., double _max_va_frac=1.)
 Construct and initialize a StrucMapper. More...
 
double cost_tol () const
 Tolerance for determining if two mapping-cost values are identical. More...
 
double xtal_tol () const
 Tolerance for initializing lattices. For now it is initialized to CASM::TOL. More...
 
double lattice_weight () const
 
double atomic_weight () const
 
void set_lattice_weight (double _lw)
 
Index lattice_transformation_range () const
 Max element considered for integer unimodular matrix transformations (which define orientation relationship of mapping) More...
 
void set_lattice_transformation_range (Index _new_range)
 Max element considered for integer unimodular matrix transformations (which define orientation relationship of mapping) More...
 
void set_symmetrize_lattice_cost (bool _sym_lat_cost)
 Flag that enables the calculation of a symmetrized lattice cost when performing the lattice maps. This cost only accounts for that part of the deformation gradient that breaks the symmetry of the parent crystal structure. More...
 
void set_symmetrize_atomic_cost (bool _sym_atomic_cost, const SymOpVector &factor_group, const std::vector< Eigen::PermutationMatrix< Eigen::Dynamic, Eigen::Dynamic, Index >> &permutation_group)
 Flag that enables the calculation of a symmetrized atomic cost when performing the atomic maps. This cost only accounts for that part of the displacement field that breaks the symmetry of the parent crystal structure. More...
 
bool symmetrize_lattice_cost () const
 
bool symmetrize_atomic_cost () const
 
double min_va_frac () const
 Returns the minimum fraction of sites allowed to be vacant in the mapping relation Vacancy fraction is used to constrain the mapping supercell search, but is only used when the supercell volume cannot is not uniquely determined by the number of each species in the child structure. More...
 
void set_min_va_frac (double _min_va)
 Sets the minimum fraction of sites allowed to be vacant in the mapping relation Vacancy fraction is used to constrain the mapping supercell search, but is only used when the supercell volume cannot is not uniquely determined by the number of each species in the child structure. More...
 
double max_va_frac () const
 Returns the maximum fraction of sites allowed to be vacant in the mapping relation Vacancy fraction is used to constrain the mapping supercell search, but is only used when the supercell volume cannot is not uniquely determined by the number of each species in the child structure. More...
 
void set_max_va_frac (double _max_va)
 Sets the maximum fraction of sites allowed to be vacant in the mapping relation Vacancy fraction is used to constrain the mapping supercell search, but is only used when the supercell volume cannot is not uniquely determined by the number of each species in the child structure. More...
 
int options () const
 returns bit flag of selected options for this StrucMapper More...
 
SimpleStructure const & parent () const
 returns reference to parent structure More...
 
void add_allowed_lattice (Lattice const &_lat)
 specify a superlattice of the parent to be searched during mapping More...
 
void clear_allowed_lattices () const
 clear the list of allowed parent superlattices; all superlattices will be generated automatically, as needed (default) More...
 
bool lattices_constrained () const
 returns true if the search of parent superlattices is constrained to a pre-specified list More...
 
void set_filter (std::function< bool(Lattice const &, Lattice const &)> _filter_f)
 specify to use filtered lattices for mapping. The filter function is of the form bool filter(parent_lattice, proposed_lattice) where parent_lattice is the primitive lattice of the parent structure, and proposed_lattice is a proposed superlattice of the parent structure More...
 
void unset_filter ()
 specify not to use filtered lattice for mapping More...
 
std::set< MappingNodemap_ideal_struc (const SimpleStructure &child_struc, Index k=1, double max_cost=StrucMapping::big_inf(), double min_cost=-TOL, bool keep_invalid=false) const
 k-best mappings of ideal child structure onto parent structure Assumes that child_struc and parent_struc have lattices related by an integer transformation so that search over lattices can be replaced with simple matrix solution (faster than typical search) Atomic positions need not be ideal More...
 
std::set< MappingNodemap_deformed_struc (const SimpleStructure &child_struc, Index k=1, double max_cost=StrucMapping::big_inf(), double min_cost=-TOL, bool keep_invalid=false, SymOpVector const &child_factor_group={SymOp::identity()}) const
 k-best mappings of arbitrary child structure onto parent structure, without simplifying assumptions If k and k+1, k+2, etc mappings have identical cost, k will be increased until k+n mapping has cost greater than mapping k More...
 
std::set< MappingNodemap_deformed_struc_impose_lattice_vols (const SimpleStructure &child_struc, Index min_vol, Index max_vol, Index k=1, double max_cost=StrucMapping::big_inf(), double min_cost=-TOL, bool keep_invalid=false, SymOpVector const &child_factor_group={SymOp::identity()}) const
 k-best mappings of arbitrary child structure onto parent structure imposes simplifying assumption that solution is a superstructure of parent structure having integer volume on range [min_vol, max_vol] If k and k+1, k+2, etc mappings have identical cost, k will be increased until k+n mapping has cost greater than mapping k More...
 
std::set< MappingNodemap_deformed_struc_impose_lattice (const SimpleStructure &child_struc, const Lattice &imposed_lat, Index k=1, double max_cost=StrucMapping::big_inf(), double min_cost=-TOL, bool keep_invalid=false, SymOpVector const &child_factor_group={SymOp::identity()}) const
 k-best mappings of arbitrary child structure onto parent structure imposes simplifying assumption that solution is THE superstructure of parent given by More...
 
std::set< MappingNodemap_deformed_struc_impose_lattice_node (const SimpleStructure &child_struc, const LatticeNode &imposed_node, Index k=1, double max_cost=StrucMapping::big_inf(), double min_cost=-TOL, bool keep_invalid=false) const
 k-best mappings of arbitrary child structure onto parent structure imposes simplifying assumption that the lattice mapping is known and specified by More...
 
Index k_best_maps_better_than (SimpleStructure const &child_struc, std::set< MappingNode > &queue, Index k=1, double max_cost=StrucMapping::big_inf(), double min_cost=-TOL, bool keep_invalid=false, bool keep_tail=false, bool no_partiton=false) const
 low-level function. Takes a queue of mappings and use it to seed a search for k-best total mappings. The seed mappings may be partial (having only LatticeNode specified, and not HungarianNode) or complete. The result is appended to More...
 
StrucMapCalculatorInterface const & calculator () const
 

Private Member Functions

std::set< MappingNode_seed_k_best_from_super_lats (SimpleStructure const &child_struc, std::vector< Lattice > const &_parent_scels, std::vector< Lattice > const &_child_scels, Index k, double max_strain_cost, double min_strain_cost, SymOpVector const &child_factor_group={SymOp::identity()}) const
 generate list of partial mapping seeds (i.e., LatticeNode only) from a list of supercells of the parent structure and a list supercells of the child structure More...
 
std::set< MappingNode_seed_from_vol_range (SimpleStructure const &child_struc, Index k, Index min_vol, Index max_vol, double max_strain_cost, double min_strain_cost, SymOpVector const &child_factor_group={SymOp::identity()}) const
 construct partial mapping nodes (with uninitialized atomic_node) based on current settings considers supercells with integer volume between min_vol and max_vol More...
 
Index _n_species (SimpleStructure const &sstruc) const
 returns number of species in a SimpleStructure given the current calculator settings. Use instead of sstruc.n_atom() for consistency More...
 
bool _filter_lat (Lattice const &_parent_lat, Lattice const &_child_lat) const
 
std::pair< Index, Index_vol_range (const SimpleStructure &child_struc) const
 
std::vector< Lattice_lattices_of_vol (Index prim_vol) const
 

Private Attributes

notstd::cloneable_ptr< StrucMapCalculatorInterfacem_calc_ptr
 
Eigen::MatrixXd m_strain_gram_mat
 
double m_lattice_weight
 
double m_max_volume_change
 
int m_options
 
double m_cost_tol
 
double m_xtal_tol
 
double m_min_va_frac
 
double m_max_va_frac
 
Index m_lattice_transformation_range
 
bool m_symmetrize_lattice_cost
 
bool m_symmetrize_atomic_cost
 
bool m_filtered
 
std::function< bool(Lattice const &, Lattice const &)> m_filter_f
 
LatMapType m_superlat_map
 Maps the supercell volume to a vector of Lattices with that volume. More...
 
LatMapType m_allowed_superlat_map
 

Member Typedef Documentation

◆ LatMapType

using CASM::xtal::StrucMapper::LatMapType = std::map<Index, std::vector<Lattice> >

Definition at line 406 of file StrucMapping.hh.

Member Enumeration Documentation

◆ Options

Enumerator
none 
strict 
robust 
sym_strain 
sym_basis 
soft_va_limit 

Definition at line 408 of file StrucMapping.hh.

Constructor & Destructor Documentation

◆ StrucMapper()

CASM::xtal::StrucMapper::StrucMapper ( StrucMapCalculatorInterface const &  _calculator,
double  _lattice_weight = 0.5,
double  _max_volume_change = 0.5,
int  _options = 0,
double  _cost_tol = TOL,
double  _min_va_frac = 0.,
double  _max_va_frac = 1. 
)

Construct and initialize a StrucMapper.

Parameters
_calculatorspecialization of base class StrucMapCalculatorInterface, which controls the way that the cost_matrix and costs are calculated, determines validity of proposed mappings, and finalizes the representation of proposed mappings
Parameters
_lattice_weightfree parameter 'w' in the cost function: total_cost = w*lattice_deformation+(1-w)*atomic_deformation can vary between 0 (completely basis-focused) and 1 (completely lattice-focused)
Parameters
_max_volume_changeconstrains the search space by assuming a limit on allowed volume change only taken into account when non-interstitial vacancies are allowed in parent structure
Parameters
_options

specify a combination of StrucMapper::Options using bitwise OR: Ex. _options=StrucMapper::robust|StrucMapper::strict Options are: 'robust': does not assume the imported structure might be ideal ('robust' is much slower for importing ideal structures, but if 'robust' is not set and a non-ideal structure is passed, this will be almost always be detected and robust methods will be used instead. Thus, 'robust' is slightly faster if imported Structures are not ideal)

'strict': prevents transformation into canonical form. Tries to preserve original orientation of imported structure if possible

Parameters
_cost_toltolerance for mapping comparisons
_min_va_fracminimum fraction of vacant sites, below this fraction a mapping will not be considered
_max_va_fracmaximum fraction of vacant sites, above this fraction a mapping will not be considered

Definition at line 399 of file StrucMapping.cc.

Member Function Documentation

◆ _filter_lat()

bool CASM::xtal::StrucMapper::_filter_lat ( Lattice const &  _parent_lat,
Lattice const &  _child_lat 
) const
inlineprivate

Definition at line 772 of file StrucMapping.hh.

◆ _lattices_of_vol()

std::vector< Lattice > CASM::xtal::StrucMapper::_lattices_of_vol ( Index  prim_vol) const
private

Definition at line 687 of file StrucMapping.cc.

◆ _n_species()

Index CASM::xtal::StrucMapper::_n_species ( SimpleStructure const &  sstruc) const
private

returns number of species in a SimpleStructure given the current calculator settings. Use instead of sstruc.n_atom() for consistency

Definition at line 427 of file StrucMapping.cc.

◆ _seed_from_vol_range()

std::set< MappingNode > CASM::xtal::StrucMapper::_seed_from_vol_range ( SimpleStructure const &  child_struc,
Index  k,
Index  min_vol,
Index  max_vol,
double  max_strain_cost,
double  min_strain_cost,
SymOpVector const &  child_factor_group = {SymOp::identity()} 
) const
private

construct partial mapping nodes (with uninitialized atomic_node) based on current settings considers supercells with integer volume between min_vol and max_vol

Definition at line 520 of file StrucMapping.cc.

◆ _seed_k_best_from_super_lats()

std::set< MappingNode > CASM::xtal::StrucMapper::_seed_k_best_from_super_lats ( SimpleStructure const &  child_struc,
std::vector< Lattice > const &  _parent_scels,
std::vector< Lattice > const &  _child_scels,
Index  k,
double  max_strain_cost,
double  min_strain_cost,
SymOpVector const &  child_factor_group = {SymOp::identity()} 
) const
private

generate list of partial mapping seeds (i.e., LatticeNode only) from a list of supercells of the parent structure and a list supercells of the child structure

Parameters
kNumber of k-best mapping relations to return
max_lattice_costSearch will terminate once no lattice mappings better than max_lattice_cost are found
min_lattice_costAll lattice mappings better than min_lattice_cost will be returned, without contributing to 'k'

Definition at line 852 of file StrucMapping.cc.

◆ _vol_range()

std::pair< Index, Index > CASM::xtal::StrucMapper::_vol_range ( const SimpleStructure child_struc) const
private

Definition at line 433 of file StrucMapping.cc.

◆ add_allowed_lattice()

void CASM::xtal::StrucMapper::add_allowed_lattice ( Lattice const &  _lat)
inline

specify a superlattice of the parent to be searched during mapping

Definition at line 557 of file StrucMapping.hh.

◆ atomic_weight()

double CASM::xtal::StrucMapper::atomic_weight ( ) const
inline

Definition at line 484 of file StrucMapping.hh.

◆ calculator()

StrucMapCalculatorInterface const& CASM::xtal::StrucMapper::calculator ( ) const
inline

Definition at line 739 of file StrucMapping.hh.

◆ clear_allowed_lattices()

void CASM::xtal::StrucMapper::clear_allowed_lattices ( ) const
inline

clear the list of allowed parent superlattices; all superlattices will be generated automatically, as needed (default)

Definition at line 565 of file StrucMapping.hh.

◆ cost_tol()

double CASM::xtal::StrucMapper::cost_tol ( ) const
inline

Tolerance for determining if two mapping-cost values are identical.

Definition at line 476 of file StrucMapping.hh.

◆ k_best_maps_better_than()

Index CASM::xtal::StrucMapper::k_best_maps_better_than ( SimpleStructure const &  child_struc,
std::set< MappingNode > &  queue,
Index  k = 1,
double  max_cost = StrucMapping::big_inf(),
double  min_cost = -TOL,
bool  keep_invalid = false,
bool  keep_tail = false,
bool  no_partiton = false 
) const

low-level function. Takes a queue of mappings and use it to seed a search for k-best total mappings. The seed mappings may be partial (having only LatticeNode specified, and not HungarianNode) or complete. The result is appended to

Parameters
queue,andpartial and invalid mappings are removed from the queue. If k and k+1, k+2, etc mappings have identical cost, k will be increased until k+n mapping has cost greater than mapping k
child_strucInput structure to be mapped onto parent structure
queuelist of partial mappings to seed search. Partial mappings are removed from list as they are searched finalized mappings are inserted into list as they are found.
Parameters
kNumber of k-best mapping relations to return
max_costSearch will terminate once no mappings better than max_cost are found
min_costAll mappings better than min_cost will be reported, without contributing to 'k'
keep_invalidIf true, invalid mappings are retained in output; otherwise they are discarded
keep_tailIf true, any partial or unresolved mappings at back of queue will be retained (they are deleted otherwise)
no_partitionIf true, search over suboptimal basis permutations will be skipped (The optimal mapping will be found, but suboptimal mappings may be excluded in the k-best case. Degenerate permutations of atoms will not be identified.)
Returns
Index specifying the number of complete, valid mappings in the solution set

Definition at line 754 of file StrucMapping.cc.

◆ lattice_transformation_range()

Index CASM::xtal::StrucMapper::lattice_transformation_range ( ) const
inline

Max element considered for integer unimodular matrix transformations (which define orientation relationship of mapping)

Definition at line 492 of file StrucMapping.hh.

◆ lattice_weight()

double CASM::xtal::StrucMapper::lattice_weight ( ) const
inline

Definition at line 482 of file StrucMapping.hh.

◆ lattices_constrained()

bool CASM::xtal::StrucMapper::lattices_constrained ( ) const
inline

returns true if the search of parent superlattices is constrained to a pre-specified list

Definition at line 569 of file StrucMapping.hh.

◆ map_deformed_struc()

std::set< MappingNode > CASM::xtal::StrucMapper::map_deformed_struc ( const SimpleStructure child_struc,
Index  k = 1,
double  max_cost = StrucMapping::big_inf(),
double  min_cost = -TOL,
bool  keep_invalid = false,
SymOpVector const &  child_factor_group = {SymOp::identity()} 
) const

k-best mappings of arbitrary child structure onto parent structure, without simplifying assumptions If k and k+1, k+2, etc mappings have identical cost, k will be increased until k+n mapping has cost greater than mapping k

Parameters
child_strucInput structure to be mapped onto parent structure
kNumber of k-best mapping relations to return
max_costSearch will terminate once no mappings better than max_cost are found
min_costAll mappings better than min_cost will be reported, without contributing to 'k'
keep_invalidIf true, invalid mappings are retained in output; otherwise they are discarded
Returns
std::set<MappingNode> a list of valid mappings, sorted first by cost, and then other attributes

Definition at line 622 of file StrucMapping.cc.

◆ map_deformed_struc_impose_lattice()

std::set< MappingNode > CASM::xtal::StrucMapper::map_deformed_struc_impose_lattice ( const SimpleStructure child_struc,
const Lattice imposed_lat,
Index  k = 1,
double  max_cost = StrucMapping::big_inf(),
double  min_cost = -TOL,
bool  keep_invalid = false,
SymOpVector const &  child_factor_group = {SymOp::identity()} 
) const

k-best mappings of arbitrary child structure onto parent structure imposes simplifying assumption that solution is THE superstructure of parent given by

Parameters
imposed_lat.Significantly faster than unconstrained search, but may miss best solution If k and k+1, k+2, etc mappings have identical cost, k will be increased until k+n mapping has cost greater than mapping k
child_strucInput structure to be mapped onto parent structure
kNumber of k-best mapping relations to return
max_costSearch will terminate once no mappings better than max_cost are found
min_costAll mappings better than min_cost will be reported, without contributing to 'k'
keep_invalidIf true, invalid mappings are retained in output; otherwise they are discarded
Returns
std::set<MappingNode> a list of valid mappings, sorted first by cost, and then other attributes

Definition at line 654 of file StrucMapping.cc.

◆ map_deformed_struc_impose_lattice_node()

std::set< MappingNode > CASM::xtal::StrucMapper::map_deformed_struc_impose_lattice_node ( const SimpleStructure child_struc,
const LatticeNode imposed_node,
Index  k = 1,
double  max_cost = StrucMapping::big_inf(),
double  min_cost = -TOL,
bool  keep_invalid = false 
) const

k-best mappings of arbitrary child structure onto parent structure imposes simplifying assumption that the lattice mapping is known and specified by

Parameters
imposed_node.This sets the lattice and setting that defines the mapping relation This assumption is much than unconstrained search, but may miss best solution If k and k+1, k+2, etc mappings have identical cost, k will be increased until k+n mapping has cost greater than mapping k
child_strucInput structure to be mapped onto parent structure
kNumber of k-best mapping relations to return
max_costSearch will terminate once no mappings better than max_cost are found
min_costAll mappings better than min_cost will be reported, without contributing to 'k'
keep_invalidIf true, invalid mappings are retained in output; otherwise they are discarded
Returns
std::set<MappingNode> a list of valid mappings, sorted first by cost, and then other attributes

Definition at line 673 of file StrucMapping.cc.

◆ map_deformed_struc_impose_lattice_vols()

std::set< MappingNode > CASM::xtal::StrucMapper::map_deformed_struc_impose_lattice_vols ( const SimpleStructure child_struc,
Index  min_vol,
Index  max_vol,
Index  k = 1,
double  max_cost = StrucMapping::big_inf(),
double  min_cost = -TOL,
bool  keep_invalid = false,
SymOpVector const &  child_factor_group = {SymOp::identity()} 
) const

k-best mappings of arbitrary child structure onto parent structure imposes simplifying assumption that solution is a superstructure of parent structure having integer volume on range [min_vol, max_vol] If k and k+1, k+2, etc mappings have identical cost, k will be increased until k+n mapping has cost greater than mapping k

Parameters
child_strucInput structure to be mapped onto parent structure
kNumber of k-best mapping relations to return
max_costSearch will terminate once no mappings better than max_cost are found
min_costAll mappings better than min_cost will be reported, without contributing to 'k'
keep_invalidIf true, invalid mappings are retained in output; otherwise they are discarded
Returns
std::set<MappingNode> a list of valid mappings, sorted first by cost, and then other attributes

Definition at line 634 of file StrucMapping.cc.

◆ map_ideal_struc()

std::set< MappingNode > CASM::xtal::StrucMapper::map_ideal_struc ( const SimpleStructure child_struc,
Index  k = 1,
double  max_cost = StrucMapping::big_inf(),
double  min_cost = -TOL,
bool  keep_invalid = false 
) const

k-best mappings of ideal child structure onto parent structure Assumes that child_struc and parent_struc have lattices related by an integer transformation so that search over lattices can be replaced with simple matrix solution (faster than typical search) Atomic positions need not be ideal

Parameters
child_strucInput structure to be mapped onto parent structure
kNumber of k-best mapping relations to return.
max_costSearch will terminate once no mappings better than max_cost are found
min_costAll mappings better than min_cost will be reported, without contributing to 'k'
Returns
std::set<MappingNode> a list of valid mappings, sorted first by cost, and then other attributes

Definition at line 573 of file StrucMapping.cc.

◆ max_va_frac()

double CASM::xtal::StrucMapper::max_va_frac ( ) const
inline

Returns the maximum fraction of sites allowed to be vacant in the mapping relation Vacancy fraction is used to constrain the mapping supercell search, but is only used when the supercell volume cannot is not uniquely determined by the number of each species in the child structure.

Definition at line 542 of file StrucMapping.hh.

◆ min_va_frac()

double CASM::xtal::StrucMapper::min_va_frac ( ) const
inline

Returns the minimum fraction of sites allowed to be vacant in the mapping relation Vacancy fraction is used to constrain the mapping supercell search, but is only used when the supercell volume cannot is not uniquely determined by the number of each species in the child structure.

Definition at line 530 of file StrucMapping.hh.

◆ options()

int CASM::xtal::StrucMapper::options ( ) const
inline

returns bit flag of selected options for this StrucMapper

Definition at line 551 of file StrucMapping.hh.

◆ parent()

SimpleStructure const & CASM::xtal::StrucMapper::parent ( ) const

returns reference to parent structure

Definition at line 514 of file StrucMapping.cc.

◆ set_filter()

void CASM::xtal::StrucMapper::set_filter ( std::function< bool(Lattice const &, Lattice const &)>  _filter_f)
inline

specify to use filtered lattices for mapping. The filter function is of the form bool filter(parent_lattice, proposed_lattice) where parent_lattice is the primitive lattice of the parent structure, and proposed_lattice is a proposed superlattice of the parent structure

Definition at line 576 of file StrucMapping.hh.

◆ set_lattice_transformation_range()

void CASM::xtal::StrucMapper::set_lattice_transformation_range ( Index  _new_range)
inline

Max element considered for integer unimodular matrix transformations (which define orientation relationship of mapping)

Definition at line 498 of file StrucMapping.hh.

◆ set_lattice_weight()

void CASM::xtal::StrucMapper::set_lattice_weight ( double  _lw)
inline

Definition at line 486 of file StrucMapping.hh.

◆ set_max_va_frac()

void CASM::xtal::StrucMapper::set_max_va_frac ( double  _max_va)
inline

Sets the maximum fraction of sites allowed to be vacant in the mapping relation Vacancy fraction is used to constrain the mapping supercell search, but is only used when the supercell volume cannot is not uniquely determined by the number of each species in the child structure.

Definition at line 548 of file StrucMapping.hh.

◆ set_min_va_frac()

void CASM::xtal::StrucMapper::set_min_va_frac ( double  _min_va)
inline

Sets the minimum fraction of sites allowed to be vacant in the mapping relation Vacancy fraction is used to constrain the mapping supercell search, but is only used when the supercell volume cannot is not uniquely determined by the number of each species in the child structure.

Definition at line 536 of file StrucMapping.hh.

◆ set_symmetrize_atomic_cost()

void CASM::xtal::StrucMapper::set_symmetrize_atomic_cost ( bool  _sym_atomic_cost,
const SymOpVector factor_group,
const std::vector< Eigen::PermutationMatrix< Eigen::Dynamic, Eigen::Dynamic, Index >> &  permutation_group 
)
inline

Flag that enables the calculation of a symmetrized atomic cost when performing the atomic maps. This cost only accounts for that part of the displacement field that breaks the symmetry of the parent crystal structure.

Definition at line 514 of file StrucMapping.hh.

◆ set_symmetrize_lattice_cost()

void CASM::xtal::StrucMapper::set_symmetrize_lattice_cost ( bool  _sym_lat_cost)
inline

Flag that enables the calculation of a symmetrized lattice cost when performing the lattice maps. This cost only accounts for that part of the deformation gradient that breaks the symmetry of the parent crystal structure.

Definition at line 506 of file StrucMapping.hh.

◆ symmetrize_atomic_cost()

bool CASM::xtal::StrucMapper::symmetrize_atomic_cost ( ) const
inline

Definition at line 524 of file StrucMapping.hh.

◆ symmetrize_lattice_cost()

bool CASM::xtal::StrucMapper::symmetrize_lattice_cost ( ) const
inline

Definition at line 523 of file StrucMapping.hh.

◆ unset_filter()

void CASM::xtal::StrucMapper::unset_filter ( )
inline

specify not to use filtered lattice for mapping

Definition at line 584 of file StrucMapping.hh.

◆ xtal_tol()

double CASM::xtal::StrucMapper::xtal_tol ( ) const
inline

Tolerance for initializing lattices. For now it is initialized to CASM::TOL.

Definition at line 480 of file StrucMapping.hh.

Member Data Documentation

◆ m_allowed_superlat_map

LatMapType CASM::xtal::StrucMapper::m_allowed_superlat_map
mutableprivate

Definition at line 801 of file StrucMapping.hh.

◆ m_calc_ptr

notstd::cloneable_ptr<StrucMapCalculatorInterface> CASM::xtal::StrucMapper::m_calc_ptr
private

Definition at line 779 of file StrucMapping.hh.

◆ m_cost_tol

double CASM::xtal::StrucMapper::m_cost_tol
private

Definition at line 786 of file StrucMapping.hh.

◆ m_filter_f

std::function<bool(Lattice const &, Lattice const &)> CASM::xtal::StrucMapper::m_filter_f
private

Definition at line 797 of file StrucMapping.hh.

◆ m_filtered

bool CASM::xtal::StrucMapper::m_filtered
private

Definition at line 796 of file StrucMapping.hh.

◆ m_lattice_transformation_range

Index CASM::xtal::StrucMapper::m_lattice_transformation_range
private

Definition at line 791 of file StrucMapping.hh.

◆ m_lattice_weight

double CASM::xtal::StrucMapper::m_lattice_weight
private

Definition at line 783 of file StrucMapping.hh.

◆ m_max_va_frac

double CASM::xtal::StrucMapper::m_max_va_frac
private

Definition at line 789 of file StrucMapping.hh.

◆ m_max_volume_change

double CASM::xtal::StrucMapper::m_max_volume_change
private

Definition at line 784 of file StrucMapping.hh.

◆ m_min_va_frac

double CASM::xtal::StrucMapper::m_min_va_frac
private

Definition at line 788 of file StrucMapping.hh.

◆ m_options

int CASM::xtal::StrucMapper::m_options
private

Definition at line 785 of file StrucMapping.hh.

◆ m_strain_gram_mat

Eigen::MatrixXd CASM::xtal::StrucMapper::m_strain_gram_mat
private

Definition at line 781 of file StrucMapping.hh.

◆ m_superlat_map

LatMapType CASM::xtal::StrucMapper::m_superlat_map
mutableprivate

Maps the supercell volume to a vector of Lattices with that volume.

Definition at line 800 of file StrucMapping.hh.

◆ m_symmetrize_atomic_cost

bool CASM::xtal::StrucMapper::m_symmetrize_atomic_cost
private

Definition at line 794 of file StrucMapping.hh.

◆ m_symmetrize_lattice_cost

bool CASM::xtal::StrucMapper::m_symmetrize_lattice_cost
private

Definition at line 793 of file StrucMapping.hh.

◆ m_xtal_tol

double CASM::xtal::StrucMapper::m_xtal_tol
private

Definition at line 787 of file StrucMapping.hh.


The documentation for this class was generated from the following files: