CASM  1.1.0
A Clusters Approach to Statistical Mechanics
Strain.cc
Go to the documentation of this file.
2 
3 #include "casm/external/Eigen/Core"
4 
5 namespace CASM {
6 namespace strain {
7 
8 template <>
9 Eigen::Matrix3d deformation_tensor_to_metric<METRIC::GREEN_LAGRANGE>(
10  const Eigen::Ref<const Eigen::Matrix3d> &deformation_tensor) {
11  const Eigen::Matrix3d &F = deformation_tensor;
12  return 0.5 * (F.transpose() * F - Eigen::MatrixXd::Identity(3, 3));
13 }
14 
15 template <>
16 Eigen::Matrix3d metric_to_deformation_tensor<METRIC::GREEN_LAGRANGE>(
17  const Eigen::Ref<const Eigen::Matrix3d> &metric_tensor) {
18  const Eigen::Matrix3d &E = metric_tensor;
19  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d> es(
20  2 * E + Eigen::MatrixXd::Identity(3, 3));
21  return es.operatorSqrt();
22 }
23 
24 //***********************************************************//
25 
26 template <>
27 Eigen::Matrix3d deformation_tensor_to_metric<METRIC::BIOT>(
28  const Eigen::Ref<const Eigen::Matrix3d> &deformation_tensor) {
29  const Eigen::Matrix3d &F = deformation_tensor;
31 }
32 
33 template <>
34 Eigen::Matrix3d metric_to_deformation_tensor<METRIC::BIOT>(
35  const Eigen::Ref<const Eigen::Matrix3d> &metric_tensor) {
36  const Eigen::Matrix3d &B = metric_tensor;
37  return B + Eigen::MatrixXd::Identity(3, 3);
38 }
39 
40 //***********************************************************//
41 
42 template <>
43 Eigen::Matrix3d deformation_tensor_to_metric<METRIC::HENCKY>(
44  const Eigen::Ref<const Eigen::Matrix3d> &deformation_tensor) {
45  const Eigen::Matrix3d &F = deformation_tensor;
46  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d> es(F.transpose() * F);
47  return es.eigenvectors() *
48  es.eigenvalues().array().log().matrix().asDiagonal() *
49  es.eigenvectors().inverse() / 2.0;
50 }
51 
52 template <>
53 Eigen::Matrix3d metric_to_deformation_tensor<METRIC::HENCKY>(
54  const Eigen::Ref<const Eigen::Matrix3d> &metric_tensor) {
55  const Eigen::Matrix3d &H = metric_tensor;
56  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d> es(H);
57  return es.eigenvectors() *
58  es.eigenvalues().array().exp().matrix().asDiagonal() *
59  es.eigenvectors().inverse();
60 }
61 
62 //***********************************************************//
63 
64 template <>
65 Eigen::Matrix3d deformation_tensor_to_metric<METRIC::EULER_ALMANSI>(
66  const Eigen::Ref<const Eigen::Matrix3d> &deformation_tensor) {
67  const Eigen::Matrix3d &F = deformation_tensor;
68  return 0.5 *
69  (Eigen::MatrixXd::Identity(3, 3) - (F * F.transpose()).inverse());
70 }
71 
72 template <>
73 Eigen::Matrix3d metric_to_deformation_tensor<METRIC::EULER_ALMANSI>(
74  const Eigen::Ref<const Eigen::Matrix3d> &metric_tensor) {
75  const Eigen::Matrix3d &A = metric_tensor;
76  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d> es(
77  Eigen::MatrixXd::Identity(3, 3) - 2 * A);
78  return es.operatorInverseSqrt();
79 }
80 
81 //***********************************************************//
82 
83 /* template <> */
84 /* Eigen::Matrix3d deformation_tensor_to_metric<METRIC::DISP_GRAD>(const
85  * Eigen::Ref<const Eigen::Matrix3d>& deformation_tensor) */
86 /* { */
87 /* const Eigen::Matrix3d& F=deformation_tensor; */
88 /* } */
89 
90 /* template <> */
91 /* Eigen::Matrix3d metric_to_deformation_tensor<METRIC::DISP_GRAD>(const
92  * Eigen::Ref<const Eigen::Matrix3d>& metric_tensor) */
93 /* { */
94 /* } */
95 
96 //***********************************************************//
97 
99  const Eigen::Ref<const Eigen::Matrix3d> &deformation_tensor) {
100  const Eigen::Matrix3d &F = deformation_tensor;
101  return F.transpose() * F;
102 }
103 
105  const Eigen::Ref<const Eigen::Matrix3d> &deformation_tensor) {
106  const Eigen::Matrix3d &F = deformation_tensor;
108  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d> eigen_solver(C);
109  return eigen_solver.operatorSqrt();
110 }
111 
112 } // namespace strain
113 } // namespace CASM
Eigen::Matrix< typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime > inverse(const Eigen::MatrixBase< Derived > &M)
Return the integer inverse matrix of an invertible integer matrix.
IdentitySymRepBuilder Identity()
Eigen::Matrix3d right_stretch_tensor(const Eigen::Ref< const Eigen::Matrix3d > &deformation_tensor)
Calculates and returns the value of U where F = R*U.
Definition: Strain.cc:104
Eigen::Matrix3d metric_tensor(const Eigen::Ref< const Eigen::Matrix3d > &deformation_tensor)
Definition: Strain.cc:98
Main CASM namespace.
Definition: APICommand.hh:8
Eigen::Matrix3d Matrix3d