CASM
AClustersApproachtoStatisticalMechanics
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules
LatticeMap.hh
Go to the documentation of this file.
1 #ifndef LATTICEMAP_HH
2 #define LATTICEMAP_HH
3 
5 
6 namespace CASM {
7 
12 
25 
26  class LatticeMap {
27  public:
28  typedef Eigen::Matrix<double, 3, 3, Eigen::DontAlign> DMatType;
29  typedef Eigen::Matrix<int, 3, 3, Eigen::DontAlign> IMatType;
30 
31  LatticeMap(const Lattice &_ideal, const Lattice &_strained, Index _num_atoms, double _tol /*= TOL*/, int _range /*= 2*/);
32  // Finds the smallest strain tensor (in terms of Frobenius norm) that deforms (*this) into a lattice symmetrically equivalent to 'strained_lattice'
33  const LatticeMap &best_strain_mapping() const;
34  const LatticeMap &next_mapping_better_than(double max_cost) const;
35  double strain_cost() const {
36  return m_cost;
37  }
38  const DMatType &matrixN() const {
39  return m_N;
40  }
41  const DMatType &matrixF() const {
42  return m_F;
43  }
44 
45  // calculated strain cost function given deformation gradient 'F' and volume of the relaxed cell
46  static double calc_strain_cost(const Eigen::Matrix3d &F, double relaxed_atomic_vol);
47  private:
48  DMatType m_L1, m_L2;
49  //Conversion matrices:
50  // m_N = m_U * m_inv_count().inverse() * m_V_inv
51  DMatType m_U, m_V_inv;
52  // m_scale = (det(m_L2)/det(m_L1))^(1/3) = det(m_F)^(1/3)
54  double m_tol;
55 
56  mutable double m_cost;
58  mutable DMatType m_F, m_N, m_cache;
59 
60  const LatticeMap &_next_mapping_better_than(double max_cost) const;
61  // use m_F and m_atomic_vol to calculate strain cost
62  double _calc_strain_cost() const;
63 
64  };
65 
67 }
68 #endif
double m_atomic_vol
Definition: LatticeMap.hh:53
Eigen::Matrix< int, 3, 3, Eigen::DontAlign > IMatType
Definition: LatticeMap.hh:29
A Counter allows looping over many incrementing variables in one loop.
Definition: Counter.hh:71
double _calc_strain_cost() const
Definition: LatticeMap.cc:126
Main CASM namespace.
Definition: complete.cpp:8
LatticeMap(const Lattice &_ideal, const Lattice &_strained, Index _num_atoms, double _tol, int _range)
Definition: LatticeMap.cc:4
double strain_cost() const
Definition: LatticeMap.hh:35
const LatticeMap & best_strain_mapping() const
Definition: LatticeMap.cc:47
EigenIndex Index
For long integer indexing:
DMatType m_cache
Definition: LatticeMap.hh:58
const DMatType & matrixF() const
Definition: LatticeMap.hh:41
const LatticeMap & _next_mapping_better_than(double max_cost) const
Definition: LatticeMap.cc:78
const DMatType & matrixN() const
Definition: LatticeMap.hh:38
Eigen::Matrix3d Matrix3d
static double calc_strain_cost(const Eigen::Matrix3d &F, double relaxed_atomic_vol)
Definition: LatticeMap.cc:118
Eigen::Matrix< double, 3, 3, Eigen::DontAlign > DMatType
Definition: LatticeMap.hh:28
const LatticeMap & next_mapping_better_than(double max_cost) const
Definition: LatticeMap.cc:70
DMatType m_V_inv
Definition: LatticeMap.hh:51
EigenCounter< IMatType > m_inv_count
Definition: LatticeMap.hh:57