CASM  1.1.0
A Clusters Approach to Statistical Mechanics
CASM::strain Namespace Reference

Enumerations

enum class  METRIC {
  GREEN_LAGRANGE , BIOT , HENCKY , EULER_ALMANSI ,
  STRETCH , DISP_GRAD
}
 

Functions

template<METRIC >
Eigen::Matrix3d deformation_tensor_to_metric (const Eigen::Ref< const Eigen::Matrix3d > &deformation_tensor)
 
template<METRIC >
Eigen::Matrix3d metric_to_deformation_tensor (const Eigen::Ref< const Eigen::Matrix3d > &metric_tensor)
 
Eigen::Matrix3d metric_tensor (const Eigen::Ref< const Eigen::Matrix3d > &deformation_tensor)
 
Eigen::Matrix3d right_stretch_tensor (const Eigen::Ref< const Eigen::Matrix3d > &deformation_tensor)
 Calculates and returns the value of U where F = R*U. More...
 
template<>
Eigen::Matrix3d deformation_tensor_to_metric< METRIC::GREEN_LAGRANGE > (const Eigen::Ref< const Eigen::Matrix3d > &deformation_tensor)
 
template<>
Eigen::Matrix3d metric_to_deformation_tensor< METRIC::GREEN_LAGRANGE > (const Eigen::Ref< const Eigen::Matrix3d > &metric_tensor)
 
template<>
Eigen::Matrix3d deformation_tensor_to_metric< METRIC::BIOT > (const Eigen::Ref< const Eigen::Matrix3d > &deformation_tensor)
 
template<>
Eigen::Matrix3d metric_to_deformation_tensor< METRIC::BIOT > (const Eigen::Ref< const Eigen::Matrix3d > &metric_tensor)
 
template<>
Eigen::Matrix3d deformation_tensor_to_metric< METRIC::HENCKY > (const Eigen::Ref< const Eigen::Matrix3d > &deformation_tensor)
 
template<>
Eigen::Matrix3d metric_to_deformation_tensor< METRIC::HENCKY > (const Eigen::Ref< const Eigen::Matrix3d > &metric_tensor)
 
template<>
Eigen::Matrix3d deformation_tensor_to_metric< METRIC::EULER_ALMANSI > (const Eigen::Ref< const Eigen::Matrix3d > &deformation_tensor)
 
template<>
Eigen::Matrix3d metric_to_deformation_tensor< METRIC::EULER_ALMANSI > (const Eigen::Ref< const Eigen::Matrix3d > &metric_tensor)
 

Enumeration Type Documentation

◆ METRIC

enum CASM::strain::METRIC
strong

Enumerates different types of available strain metrics, which can be used to specify conversion between deformation tensors and the desired metric:

F: deformation tensor (F = R * U) R: rotation tensor U: stretch tensor C: F^{T} * F

GREEN_LAGRANGE = 1/2 * (F^{T} F - I) BIOT = (U-I) HENCKY = log(C)/2 EULER_ALMANSI = (I-(F F^{T})^(-1))/2

Enumerator
GREEN_LAGRANGE 
BIOT 
HENCKY 
EULER_ALMANSI 
STRETCH 
DISP_GRAD 

Definition at line 24 of file Strain.hh.

Function Documentation

◆ deformation_tensor_to_metric()

template<METRIC >
Eigen::Matrix3d CASM::strain::deformation_tensor_to_metric ( const Eigen::Ref< const Eigen::Matrix3d > &  deformation_tensor)

Convert a deformation tensor to any metric tensor, specified by the template parameter. For example, in order to turn a deformation tensor F into a Green Lagrange metric tensor GL: Eigen::Matrix3d GL=deformation_tensor_to_metric<METRIC::GREEN_LAGRANGE>(F);

◆ deformation_tensor_to_metric< METRIC::BIOT >()

template<>
Eigen::Matrix3d CASM::strain::deformation_tensor_to_metric< METRIC::BIOT > ( const Eigen::Ref< const Eigen::Matrix3d > &  deformation_tensor)

Definition at line 27 of file Strain.cc.

◆ deformation_tensor_to_metric< METRIC::EULER_ALMANSI >()

template<>
Eigen::Matrix3d CASM::strain::deformation_tensor_to_metric< METRIC::EULER_ALMANSI > ( const Eigen::Ref< const Eigen::Matrix3d > &  deformation_tensor)

Definition at line 65 of file Strain.cc.

◆ deformation_tensor_to_metric< METRIC::GREEN_LAGRANGE >()

template<>
Eigen::Matrix3d CASM::strain::deformation_tensor_to_metric< METRIC::GREEN_LAGRANGE > ( const Eigen::Ref< const Eigen::Matrix3d > &  deformation_tensor)

Definition at line 9 of file Strain.cc.

◆ deformation_tensor_to_metric< METRIC::HENCKY >()

template<>
Eigen::Matrix3d CASM::strain::deformation_tensor_to_metric< METRIC::HENCKY > ( const Eigen::Ref< const Eigen::Matrix3d > &  deformation_tensor)

Definition at line 43 of file Strain.cc.

◆ metric_tensor()

Eigen::Matrix3d CASM::strain::metric_tensor ( const Eigen::Ref< const Eigen::Matrix3d > &  deformation_tensor)

Calculates the metric tensor of the the deformation gradient as F^{T}F

Definition at line 98 of file Strain.cc.

◆ metric_to_deformation_tensor()

template<METRIC >
Eigen::Matrix3d CASM::strain::metric_to_deformation_tensor ( const Eigen::Ref< const Eigen::Matrix3d > &  metric_tensor)

Convert a strain metric tensor to the deformation tensor, where the given strain metric is specified with the template parameter. For example, to convert a Green Lagrange strain tensor GL to the deformation tensor F: Eigen::Matrix3d F=metric_to_deformation_tensor<METRIC::GREEN_LAGRANGE>(GL);

◆ metric_to_deformation_tensor< METRIC::BIOT >()

template<>
Eigen::Matrix3d CASM::strain::metric_to_deformation_tensor< METRIC::BIOT > ( const Eigen::Ref< const Eigen::Matrix3d > &  metric_tensor)

Definition at line 34 of file Strain.cc.

◆ metric_to_deformation_tensor< METRIC::EULER_ALMANSI >()

template<>
Eigen::Matrix3d CASM::strain::metric_to_deformation_tensor< METRIC::EULER_ALMANSI > ( const Eigen::Ref< const Eigen::Matrix3d > &  metric_tensor)

Definition at line 73 of file Strain.cc.

◆ metric_to_deformation_tensor< METRIC::GREEN_LAGRANGE >()

template<>
Eigen::Matrix3d CASM::strain::metric_to_deformation_tensor< METRIC::GREEN_LAGRANGE > ( const Eigen::Ref< const Eigen::Matrix3d > &  metric_tensor)

Definition at line 16 of file Strain.cc.

◆ metric_to_deformation_tensor< METRIC::HENCKY >()

template<>
Eigen::Matrix3d CASM::strain::metric_to_deformation_tensor< METRIC::HENCKY > ( const Eigen::Ref< const Eigen::Matrix3d > &  metric_tensor)

Definition at line 53 of file Strain.cc.

◆ right_stretch_tensor()

Eigen::Matrix3d CASM::strain::right_stretch_tensor ( const Eigen::Ref< const Eigen::Matrix3d > &  deformation_tensor)

Calculates and returns the value of U where F = R*U.

Definition at line 104 of file Strain.cc.