CASM
AClustersApproachtoStatisticalMechanics
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules
StrainConverter.hh
Go to the documentation of this file.
1 #ifndef STRAINCONVERTER_HH
2 #define STRAINCONVERTER_HH
3 
6 
7 
8 //--------------------------------------------------
9 // STRAIN CLASS
10 // NAMING CONVENTION:
11 // - F: deformation tensor (F = R * U)
12 // - E: strain metric (calculated based on what STRAIN_METRIC_MODE is)
13 // - R: rotation tensor
14 // - U: stretch tensor
15 // - C: F^{T} * F
16 // - unrolled_E: the 6 unique components of E ordered in an order
17 // specified by order_strain
18 // - sop: strain order parameter (sop = sop_transf_mat * unrolled_E)
19 
20 namespace CASM {
21  class SymGroup;
22 
23  typedef Eigen::VectorXd VectorXd;
26 
27  //Enum type that will define how the strain metrics used to
28  //calculate strain order parameters are going to be calculated
30 
31 
33  public:
34  // Static methods for conversion from F to various strain metrics
35 
37  static Matrix3d green_lagrange(const Matrix3d &F);
38  static Matrix3d green_lagrange_to_F(const Matrix3d &E);
39 
41  static Matrix3d biot(const Matrix3d &F);
42  static Matrix3d biot_to_F(const Matrix3d &B);
43 
45  static Matrix3d hencky(const Matrix3d &F);
46  static Matrix3d hencky_to_F(const Matrix3d &H);
47 
49  static Matrix3d euler_almansi(const Matrix3d &F);
50  static Matrix3d euler_almansi_to_F(const Matrix3d &A);
51 
53  static Matrix3d disp_grad(const Matrix3d &F);
54  static Matrix3d disp_grad_to_F(const Matrix3d &F) {
55  return F;
56  }
57 
58  //-------------------------------------------------
59  // Routines that calculate derived quantities given the
60  // deformation tensor
61  static Matrix3d metric_tensor(const Matrix3d &F);
62  static Matrix3d right_stretch_tensor(Matrix3d &C, Eigen::Vector3d &eigenvalues,
63  Matrix3d &eigenvectors, const Matrix3d &F);
64  static Matrix3d right_stretch_tensor(const Matrix3d &F);
65 
66  static Matrix3d strain_metric(const Matrix3d &F, STRAIN_METRIC MODE);
67  //-------------------------------------------------
68 
69  StrainConverter(bool override_flag = false) {
70  if(!override_flag) {
71  std::cerr << "WARNING in CASM::StrainConverter you are calling the default constructor" << std::endl;
72  std::cerr << "This is going to initialize a \"default\" sop_transf_mat and order_strain" << std::endl;
73  std::cerr << "The Green-Lagrange strain metric will be used for the calculation of all deformation metrics" << std::endl;
74  std::cerr << "PLEASE USE WITH CAUTION" << std::endl;
75  }
80  };
81 
82  StrainConverter(const STRAIN_METRIC &_MODE, const MatrixXd &_sop_transf_mat,
83  const Array< Array<int> > &_order_strain) :
84  STRAIN_METRIC_MODE(_MODE), m_sop_transf_mat(_sop_transf_mat),
85  m_order_strain(_order_strain) {
86  if(_MODE == GREEN_LAGRANGE)
88  else if(_MODE == BIOT)
90  else if(_MODE == HENCKY)
92  else if(_MODE == EULER_ALMANSI)
94  else if(_MODE == DISP_GRAD)
96  }
97 
98  StrainConverter(const std::string &mode_name) :
99  StrainConverter(true) {
100  set_mode(mode_name);
101  }
102 
103  Index dim() const {
104  return m_order_strain.size();
105  }
106 
107  //-------------------------------------------------
109  Matrix3d strain_metric(const Matrix3d &F) const;
110  Matrix3d strain_metric_to_F(const Matrix3d &E) const;
111 
113  VectorXd unroll_E(const Matrix3d &E) const;
114  Matrix3d rollup_E(const VectorXd &_unrolled_E) const;
115 
116  VectorXd unrolled_strain_metric(const Matrix3d &F) const;
117  Matrix3d unrolled_strain_metric_to_F(const VectorXd &E) const;
118 
119  VectorXd sop(Matrix3d &E, Matrix3d &C, Matrix3d &U, Eigen::Vector3d &eigenvalues,
120  Matrix3d &eigenvectors, const Matrix3d &F) const;
121  VectorXd sop(Matrix3d &E, Matrix3d &C, Matrix3d &U, Eigen::Vector3d &eigenvalues,
122  Matrix3d &eigenvectors, const Matrix3d &F, STRAIN_METRIC MODE) const;
123 
124  //--------------------------------------------------
125  // Routines that set the internal parameters of the
126  // class
127  void set_mode(const std::string &mode_name);
128 
130  return m_sop_transf_mat;
131  }
132 
134  return m_symrep_ID;
135  }
136 
137  std::vector<Eigen::MatrixXd> irreducible_sop_wedges(const SymGroup &pg, std::vector<Index> &multiplicities);
138 
139  std::vector<Eigen::MatrixXd> irreducible_wedges(const SymGroup &pg, std::vector<Index> &multiplicities);
140 
141  void set_symmetrized_sop(const SymGroup &pg);
142 
146 
147  // Need to rewrite these routines to write out only the settings
148  // jsonParser &to_json(const StrainConverter &strain, jsonParser &json);
149  // void from_json(StrainConverter &strain, const jsonParser &json);
150 
151  private:
152  STRAIN_METRIC STRAIN_METRIC_MODE; //set the mode when you initialize your PRIM
153  MatrixXd m_sop_transf_mat; //Use as sop = sop_transf_mat * unrolled_E
154  // unrolled_E[i]= m_weight_strain[i]*strain(m_order_strain[i][0], m_order_strain[i][1])
155  Array< Array<int> > m_order_strain; //lists the order of strains to list unrolled_E
156  Array<double> m_weight_strain; // weights for the elements of unrolled_E
157 
159  // typedef MetricFuncPtr for method pointers that take displacement gradient tensor as argument
160  typedef Matrix3d(*MetricFuncPtr)(const Matrix3d &);
163 
164  };
165 }
166 #endif
Eigen::MatrixXd MatrixXd
void set_mode(const std::string &mode_name)
const Eigen::MatrixXd & sop_transf_mat() const
static Matrix3d hencky(const Matrix3d &F)
HENCKY = log(C)/2.
static Matrix3d biot(const Matrix3d &F)
BIOT = (U-I)
static Matrix3d green_lagrange(const Matrix3d &F)
GREEN_LAGRANGE = (C-I)/2.
MetricFuncPtr curr_inv_metric_func
StrainConverter(const STRAIN_METRIC &_MODE, const MatrixXd &_sop_transf_mat, const Array< Array< int > > &_order_strain)
Type-safe ID object for communicating and accessing Symmetry representation info. ...
static Matrix3d disp_grad_to_F(const Matrix3d &F)
static Matrix3d strain_metric(const Matrix3d &F, STRAIN_METRIC MODE)
static Matrix3d hencky_to_F(const Matrix3d &H)
HENCKY = log(C)/2.
StrainConverter(const std::string &mode_name)
Matrix3d rollup_E(const VectorXd &_unrolled_E) const
Matrix3d unrolled_strain_metric_to_F(const VectorXd &E) const
VectorXd unroll_E(const Matrix3d &E) const
Unrolls the green-lagrange metric ( or any symmetric metric)
Main CASM namespace.
Definition: complete.cpp:8
Matrix3d strain_metric_to_F(const Matrix3d &E) const
SymGroupRepID symrep_ID() const
Array< Array< int > > m_order_strain
static Matrix3d euler_almansi_to_F(const Matrix3d &A)
EULER_ALMANSI = (I-(F F^{T})^(-1))/2.
SymGroup is a collection of symmetry operations that satisfy the group property The symmetry operatio...
Definition: SymGroup.hh:33
static Matrix3d euler_almansi(const Matrix3d &F)
EULER_ALMANSI = (I-(F F^{T})^(-1))/2.
static Matrix3d disp_grad(const Matrix3d &F)
DISP_GRAD = F.
EigenIndex Index
For long integer indexing:
STRAIN_METRIC STRAIN_METRIC_MODE
Eigen::VectorXd VectorXd
void set_symmetrized_sop(const SymGroup &pg)
Array< double > m_weight_strain
VectorXd sop(Matrix3d &E, Matrix3d &C, Matrix3d &U, Eigen::Vector3d &eigenvalues, Matrix3d &eigenvectors, const Matrix3d &F) const
StrainConverter(bool override_flag=false)
Eigen::Matrix3d Matrix3d
std::vector< Eigen::MatrixXd > irreducible_wedges(const SymGroup &pg, std::vector< Index > &multiplicities)
void set_conventional_order_unsymmetric()
static Matrix3d metric_tensor(const Matrix3d &F)
static Matrix3d green_lagrange_to_F(const Matrix3d &E)
GREEN_LAGRANGE = 1/2 * (F^{T} F - I)
static Matrix3d biot_to_F(const Matrix3d &B)
BIOT = (U-I)
VectorXd unrolled_strain_metric(const Matrix3d &F) const
static Matrix3d right_stretch_tensor(Matrix3d &C, Eigen::Vector3d &eigenvalues, Matrix3d &eigenvectors, const Matrix3d &F)
std::vector< Eigen::MatrixXd > irreducible_sop_wedges(const SymGroup &pg, std::vector< Index > &multiplicities)
Basic std::vector like container (deprecated)
Matrix3d(* MetricFuncPtr)(const Matrix3d &)
MetricFuncPtr curr_metric_func