CASM
1.1.0
A Clusters Approach to Statistical Mechanics
|
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) |
|
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 |
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);
Eigen::Matrix3d CASM::strain::deformation_tensor_to_metric< METRIC::BIOT > | ( | const Eigen::Ref< const Eigen::Matrix3d > & | deformation_tensor | ) |
Eigen::Matrix3d CASM::strain::deformation_tensor_to_metric< METRIC::EULER_ALMANSI > | ( | const Eigen::Ref< const Eigen::Matrix3d > & | deformation_tensor | ) |
Eigen::Matrix3d CASM::strain::deformation_tensor_to_metric< METRIC::GREEN_LAGRANGE > | ( | const Eigen::Ref< const Eigen::Matrix3d > & | deformation_tensor | ) |
Eigen::Matrix3d CASM::strain::deformation_tensor_to_metric< METRIC::HENCKY > | ( | const Eigen::Ref< const Eigen::Matrix3d > & | deformation_tensor | ) |
Eigen::Matrix3d CASM::strain::metric_tensor | ( | const Eigen::Ref< const Eigen::Matrix3d > & | deformation_tensor | ) |
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);
Eigen::Matrix3d CASM::strain::metric_to_deformation_tensor< METRIC::BIOT > | ( | const Eigen::Ref< const Eigen::Matrix3d > & | metric_tensor | ) |
Eigen::Matrix3d CASM::strain::metric_to_deformation_tensor< METRIC::EULER_ALMANSI > | ( | const Eigen::Ref< const Eigen::Matrix3d > & | metric_tensor | ) |
Eigen::Matrix3d CASM::strain::metric_to_deformation_tensor< METRIC::GREEN_LAGRANGE > | ( | const Eigen::Ref< const Eigen::Matrix3d > & | metric_tensor | ) |
Eigen::Matrix3d CASM::strain::metric_to_deformation_tensor< METRIC::HENCKY > | ( | const Eigen::Ref< const Eigen::Matrix3d > & | metric_tensor | ) |