CASM  1.1.0
A Clusters Approach to Statistical Mechanics
LatticeMap.hh
Go to the documentation of this file.
1 #ifndef LATTICEMAP_HH
2 #define LATTICEMAP_HH
3 
7 
8 namespace CASM {
9 namespace xtal {
10 
11 struct SymOp;
12 typedef std::vector<SymOp> SymOpVector;
13 
19  public:
20  StrainCostCalculator(Eigen::Ref<const Eigen::MatrixXd> const
21  &strain_gram_mat = Eigen::MatrixXd::Identity(9, 9));
22 
23  //\brief Isotropic strain cost, without gram matrix
24  static double isotropic_strain_cost(
25  Eigen::Matrix3d const &_deformation_gradient);
26 
27  //\brief Isotropic strain cost, without gram matrix
28  static double isotropic_strain_cost(
29  Eigen::Matrix3d const &_deformation_gradient, double _vol_factor);
30 
31  // \brief Volumetric factor :
32  // pow(abs(_deformation_gradient.determinant()),1./3.), used to normalize the
33  // strain cost to make it volume-independent
34  static double vol_factor(Eigen::Matrix3d const &_deformation_gradient) {
35  return pow(std::abs(_deformation_gradient.determinant()), 1. / 3.);
36  }
37 
38  //\brief Anisotropic strain cost; utilizes stored gram matrix to compute
39  // strain cost
40  double strain_cost(Eigen::Matrix3d const &_deformation_gradient) const;
41 
42  //\brief Anisotropic strain cost; utilizes stored gram matrix to compute
43  // strain cost
44  double strain_cost(Eigen::Matrix3d const &_deformation_gradient,
45  double _vol_factor) const;
46 
47  //\brief Symmetrized strain cost; Utilizes the parent point group symmetry to
48  // calculate only the symmetry breaking lattice cost
49  double strain_cost(Eigen::Matrix3d const &_deformation_gradient,
50  SymOpVector const &parent_point_group) const;
51 
52  private:
54  bool m_sym_cost;
55 
58 };
59 
77 
78 class LatticeMap {
79  public:
80  typedef Eigen::Matrix<double, 3, 3, Eigen::DontAlign> DMatType;
81  typedef Eigen::Matrix<int, 3, 3, Eigen::DontAlign> IMatType;
82 
83  LatticeMap(Lattice const &_parent, Lattice const &_child, Index _num_atoms,
84  int _range, SymOpVector const &_parent_point_group,
85  SymOpVector const &_child_point_group,
86  Eigen::Ref<const Eigen::MatrixXd> const &strain_gram_mat =
88  double _init_better_than = 1e20,
89  bool _symmetrize_strain_cost = false, double _xtal_tol = TOL);
90 
91  LatticeMap(Eigen::Ref<const DMatType> const &_parent,
92  Eigen::Ref<const DMatType> const &_child, Index _num_atoms,
93  int _range, SymOpVector const &_parent_point_group,
94  SymOpVector const &_child_point_group,
95  Eigen::Ref<const Eigen::MatrixXd> const &strain_gram_mat =
97  double _init_better_than = 1e20,
98  bool _symmetrize_strain_cost = false, double _xtal_tol = TOL);
99 
100  void reset(double _better_than = 1e20);
101 
102  // Finds the smallest strain tensor (in terms of Frobenius norm) that deforms
103  // (*this) into a lattice symmetrically equivalent to 'child_lattice'
104  LatticeMap const &best_strain_mapping() const;
105 
106  LatticeMap const &next_mapping_better_than(double max_cost) const;
107 
108  double strain_cost() const { return m_cost; }
109 
111 
112  const DMatType &matrixN() const { return m_N; }
113 
115  return m_deformation_gradient;
116  }
117 
118  const DMatType &parent_matrix() const { return m_parent; }
119 
120  const DMatType &child_matrix() const { return m_child; }
121 
122  double xtal_tol() const { return m_xtal_tol; }
124 
125  private:
127  // Conversion matrices:
128  // m_N = m_U * m_inv_count().inverse() * m_V_inv
130 
132 
133  // m_scale = (det(m_child)/det(m_parent))^(2/3) =
134  // det(m_deformation_gradient)^(2/3)
135  double m_vol_factor;
136 
137  int m_range;
138 
139  // pointer to static list of unimodular matrices
140  std::vector<Eigen::Matrix3i> const *m_mvec_ptr;
141 
142  // parent point group matrices, in fractional coordinates
143  std::vector<Eigen::Matrix3i> m_parent_fsym_mats;
144 
145  // parent point group in cartesian coordinates
147 
148  // child point group matrices, in fractional coordinates
149  std::vector<Eigen::Matrix3i> m_child_fsym_mats;
150 
151  // flag indicating if the symmetrized strain cost should be used while
152  // searching for the best lattice maps
154  double m_xtal_tol;
155 
156  mutable double m_cost;
157  mutable Index m_currmat;
160 
163  // We treat the unimodular matrices as the inverse of the transformation
164  // matrices that we are actively considering, allowing fewer matrix inversions
165  Eigen::Matrix3i const &inv_mat() const { return (*m_mvec_ptr)[m_currmat]; }
166 
168  Index n_mat() const { return m_mvec_ptr->size(); }
169 
171  bool _check_canonical() const;
172 
173  LatticeMap const &_next_mapping_better_than(double max_cost) const;
174 
175  // use m_deformation_gradient to calculate strain cost
176  double _calc_strain_cost() const;
177 };
178 
180 } // namespace xtal
181 } // namespace CASM
182 #endif
Eigen::Matrix< double, 3, 3, Eigen::DontAlign > DMatType
Definition: LatticeMap.hh:80
SymOpVector m_parent_point_group
Definition: LatticeMap.hh:146
Eigen::Matrix< int, 3, 3, Eigen::DontAlign > IMatType
Definition: LatticeMap.hh:81
double strain_cost() const
Definition: LatticeMap.hh:108
LatticeMap(Lattice const &_parent, Lattice const &_child, Index _num_atoms, int _range, SymOpVector const &_parent_point_group, SymOpVector const &_child_point_group, Eigen::Ref< const Eigen::MatrixXd > const &strain_gram_mat=Eigen::MatrixXd::Identity(9, 9), double _init_better_than=1e20, bool _symmetrize_strain_cost=false, double _xtal_tol=TOL)
Definition: LatticeMap.cc:129
Eigen::Matrix3i const & inv_mat() const
Returns the inverse of the current transformation matrix under consideration.
Definition: LatticeMap.hh:165
double xtal_tol() const
Definition: LatticeMap.hh:122
double _calc_strain_cost() const
const DMatType & parent_matrix() const
Definition: LatticeMap.hh:118
const DMatType & deformation_gradient() const
Definition: LatticeMap.hh:114
StrainCostCalculator m_calc
Definition: LatticeMap.hh:131
LatticeMap(Eigen::Ref< const DMatType > const &_parent, Eigen::Ref< const DMatType > const &_child, Index _num_atoms, int _range, SymOpVector const &_parent_point_group, SymOpVector const &_child_point_group, Eigen::Ref< const Eigen::MatrixXd > const &strain_gram_mat=Eigen::MatrixXd::Identity(9, 9), double _init_better_than=1e20, bool _symmetrize_strain_cost=false, double _xtal_tol=TOL)
LatticeMap const & next_mapping_better_than(double max_cost) const
Definition: LatticeMap.cc:306
bool symmetrize_strain_cost() const
Definition: LatticeMap.hh:123
Index n_mat() const
Number of unimodular matrices.
Definition: LatticeMap.hh:168
double calc_strain_cost(const Eigen::Matrix3d &deformation_gradient) const
Definition: LatticeMap.cc:296
void reset(double _better_than=1e20)
Definition: LatticeMap.cc:219
LatticeMap const & best_strain_mapping() const
Definition: LatticeMap.cc:275
DMatType m_deformation_gradient
Definition: LatticeMap.hh:158
const DMatType & child_matrix() const
Definition: LatticeMap.hh:120
const DMatType & matrixN() const
Definition: LatticeMap.hh:112
std::vector< Eigen::Matrix3i > const * m_mvec_ptr
Definition: LatticeMap.hh:140
std::vector< Eigen::Matrix3i > m_child_fsym_mats
Definition: LatticeMap.hh:149
LatticeMap const & _next_mapping_better_than(double max_cost) const
Definition: LatticeMap.cc:315
std::vector< Eigen::Matrix3i > m_parent_fsym_mats
Definition: LatticeMap.hh:143
bool _check_canonical() const
Returns true if current transformation is the canonical equivalent.
Definition: LatticeMap.cc:359
static double isotropic_strain_cost(Eigen::Matrix3d const &_deformation_gradient)
Definition: LatticeMap.cc:65
StrainCostCalculator(Eigen::Ref< const Eigen::MatrixXd > const &strain_gram_mat=Eigen::MatrixXd::Identity(9, 9))
Definition: LatticeMap.cc:14
double strain_cost(Eigen::Matrix3d const &_deformation_gradient) const
Definition: LatticeMap.cc:107
static double vol_factor(Eigen::Matrix3d const &_deformation_gradient)
Definition: LatticeMap.hh:34
IdentitySymRepBuilder Identity()
std::vector< SymOp > SymOpVector
Main CASM namespace.
Definition: APICommand.hh:8
Eigen::MatrixXd MatrixXd
Eigen::Matrix3d Matrix3d
const double TOL
Definition: definitions.hh:30
INDEX_TYPE Index
For long integer indexing:
Definition: definitions.hh:39