![]() |
FrontISTR
5.7.0
Large-scale structural analysis program with finit element method
|
This module provides functions for hyperelastic calculation. More...
Functions/Subroutines | |
subroutine | cderiv (matl, sectType, ctn, itn, inv1b, inv2b, inv3b, dibdc, d2ibdc2, strain, direction, inv4b, dibdc_ani, d2ibdc2_ani) |
This subroutine calculates derivative of the invariant with respect to Cauchy-Green tensor. More... | |
subroutine | cderiv_aniso (ctn, inv3, didc_3, d2idc2_3, direction, inv4b, dibdc_ani, d2ibdc2_ani) |
This subroutine calculates derivative of the 4th reduced invariant I4 with respect to Cauchy-Green tensor where I4 = I_3^{-1/3}*a_i*C_{ij}*a_j and a_i = direc(i). More... | |
subroutine | calelasticarrudaboyce (matl, sectType, cijkl, strain) |
This subroutine provides elastic tangent coefficient for Arruda-Boyce hyperelastic material. More... | |
subroutine | calupdateelasticarrudaboyce (matl, sectType, dstrain, dstress) |
This subroutine provides to update stress and strain for Arrude-Royce material. More... | |
subroutine | calelasticmooneyrivlin (matl, sectType, cijkl, strain, hdflag) |
This subroutine provides elastic tangent coefficient for Mooney-Rivlin hyperelastic material. More... | |
subroutine | calupdateelasticmooneyrivlin (matl, sectType, strain, stress, strain_energy, hdflag) |
This subroutine provides to update stress and strain for Mooney-Rivlin material. More... | |
subroutine | calelasticmooneyrivlinaniso (matl, sectType, cijkl, strain, cdsys, hdflag) |
This subroutine provides elastic tangent coefficient for anisotropic Mooney-Rivlin hyperelastic material. More... | |
subroutine | calupdateelasticmooneyrivlinaniso (matl, sectType, strain, stress, cdsys, hdflag) |
This subroutine provides to update stress and strain for anisotropic Mooney-Rivlin material. More... | |
This module provides functions for hyperelastic calculation.
subroutine mhyperelastic::calelasticarrudaboyce | ( | type( tmaterial ), intent(in) | matl, |
integer, intent(in) | sectType, | ||
real(kind=kreal), dimension(3,3,3,3), intent(out) | cijkl, | ||
real(kind=kreal), dimension(6), intent(in) | strain | ||
) |
This subroutine provides elastic tangent coefficient for Arruda-Boyce hyperelastic material.
[in] | matl | material rpoperties |
[in] | secttype | not used currently |
[out] | cijkl | constitutive relation |
[in] | strain | Cauchy-Lagrange strain tensor |
Definition at line 198 of file Hyperelastic.f90.
subroutine mhyperelastic::calelasticmooneyrivlin | ( | type( tmaterial ), intent(in) | matl, |
integer, intent(in) | sectType, | ||
real(kind=kreal), dimension(3,3,3,3), intent(out) | cijkl, | ||
real(kind=kreal), dimension(6), intent(in) | strain, | ||
integer(kind=kint), intent(in), optional | hdflag | ||
) |
This subroutine provides elastic tangent coefficient for Mooney-Rivlin hyperelastic material.
[in] | matl | material rpoperties |
[in] | secttype | not used curr |
[out] | cijkl | constitutive relation |
[in] | strain | Cauchy-Lagrange strain tensor |
Definition at line 287 of file Hyperelastic.f90.
subroutine mhyperelastic::calelasticmooneyrivlinaniso | ( | type( tmaterial ), intent(in) | matl, |
integer, intent(in) | sectType, | ||
real(kind=kreal), dimension(3,3,3,3), intent(out) | cijkl, | ||
real(kind=kreal), dimension(6), intent(in) | strain, | ||
real(kind=kreal), dimension(3,3), intent(in) | cdsys, | ||
integer(kind=kint), intent(in), optional | hdflag | ||
) |
This subroutine provides elastic tangent coefficient for anisotropic Mooney-Rivlin hyperelastic material.
[in] | matl | material rpoperties |
[in] | secttype | not used curr |
[out] | cijkl | constitutive relation |
[in] | strain | Cauchy-Lagrange strain tensor |
Definition at line 387 of file Hyperelastic.f90.
subroutine mhyperelastic::calupdateelasticarrudaboyce | ( | type( tmaterial ), intent(in) | matl, |
integer, intent(in) | sectType, | ||
real(kind=kreal), dimension(6), intent(in) | dstrain, | ||
real(kind=kreal), dimension(6), intent(out) | dstress | ||
) |
This subroutine provides to update stress and strain for Arrude-Royce material.
Definition at line 242 of file Hyperelastic.f90.
subroutine mhyperelastic::calupdateelasticmooneyrivlin | ( | type( tmaterial ), intent(in) | matl, |
integer, intent(in) | sectType, | ||
real(kind=kreal), dimension(6), intent(in) | strain, | ||
real(kind=kreal), dimension(6), intent(out) | stress, | ||
real(kind=kreal), intent(out) | strain_energy, | ||
integer(kind=kint), intent(in), optional | hdflag | ||
) |
This subroutine provides to update stress and strain for Mooney-Rivlin material.
[in] | matl | material properties |
[in] | secttype | not used currently |
[out] | stress | 2nd Piola-Kirchhoff stress |
[in] | strain | Green-Lagrangen strain |
[out] | strain_energy | strain energy |
Definition at line 334 of file Hyperelastic.f90.
subroutine mhyperelastic::calupdateelasticmooneyrivlinaniso | ( | type( tmaterial ), intent(in) | matl, |
integer, intent(in) | sectType, | ||
real(kind=kreal), dimension(6), intent(in) | strain, | ||
real(kind=kreal), dimension(6), intent(out) | stress, | ||
real(kind=kreal), dimension(3,3), intent(in) | cdsys, | ||
integer(kind=kint), intent(in), optional | hdflag | ||
) |
This subroutine provides to update stress and strain for anisotropic Mooney-Rivlin material.
[in] | matl | material properties |
[in] | secttype | not used currently |
[out] | stress | 2nd Piola-Kirchhoff stress |
[in] | strain | Green-Lagrangen strain |
Definition at line 440 of file Hyperelastic.f90.
subroutine mhyperelastic::cderiv | ( | type( tmaterial ), intent(in) | matl, |
integer, intent(in) | sectType, | ||
real(kind=kreal), dimension(3,3), intent(out) | ctn, | ||
real(kind=kreal), dimension(3,3), intent(out) | itn, | ||
real(kind=kreal), intent(out) | inv1b, | ||
real(kind=kreal), intent(out) | inv2b, | ||
real(kind=kreal), intent(out) | inv3b, | ||
real(kind=kreal), dimension(3,3,3), intent(out) | dibdc, | ||
real(kind=kreal), dimension(3,3,3,3,3), intent(out) | d2ibdc2, | ||
real(kind=kreal), dimension(6), intent(in) | strain, | ||
real(kind=kreal), dimension(3), intent(in), optional | direction, | ||
real(kind=kreal), intent(out), optional | inv4b, | ||
real(kind=kreal), dimension(3,3), intent(out), optional | dibdc_ani, | ||
real(kind=kreal), dimension(3,3,3,3), intent(out), optional | d2ibdc2_ani | ||
) |
This subroutine calculates derivative of the invariant with respect to Cauchy-Green tensor.
[in] | matl | material rpoperties |
[in] | secttype | not used currently |
[out] | inv1b | invariants |
[out] | inv2b | invariants |
[out] | inv3b | invariants |
[out] | dibdc | derivative of the invariant with respect to c(i,j) |
[out] | d2ibdc2 | derivative of the invariant with respect to c(i,j) |
[in] | strain | Cauchy-Lagrange strain tensor |
[out] | ctn | right Cauchy-Green deformation tensor |
[out] | itn | identity tensor |
[in] | direction | direction vector of anisotropic materials |
[out] | inv4b | invariants |
[out] | dibdc_ani | derivative of the 4th anisotropic invariant |
[out] | d2ibdc2_ani | 2nd derivative of the 4th anisotropic invariant |
Definition at line 16 of file Hyperelastic.f90.
subroutine mhyperelastic::cderiv_aniso | ( | real(kind=kreal), dimension(3,3), intent(in) | ctn, |
real(kind=kreal), intent(in) | inv3, | ||
real(kind=kreal), dimension(3,3), intent(in) | didc_3, | ||
real(kind=kreal), dimension(3,3,3,3), intent(in) | d2idc2_3, | ||
real(kind=kreal), dimension(3), intent(in) | direction, | ||
real(kind=kreal), intent(out) | inv4b, | ||
real(kind=kreal), dimension(3,3), intent(out) | dibdc_ani, | ||
real(kind=kreal), dimension(3,3,3,3), intent(out) | d2ibdc2_ani | ||
) |
This subroutine calculates derivative of the 4th reduced invariant I4 with respect to Cauchy-Green tensor where I4 = I_3^{-1/3}*a_i*C_{ij}*a_j and a_i = direc(i).
[in] | ctn | right Cauchy-Green deformation tensor |
[in] | inv3 | 3rd invariants |
[in] | didc_3 | derivative of the 3rd invariant with respect to c(i,j) |
[in] | d2idc2_3 | derivative of the 3rd invariant with respect to c(i,j) |
[in] | direction | direction vector of anisotropic materials |
[out] | inv4b | 4th reduced invariant I4 |
[out] | dibdc_ani | derivative of the 4th anisotropic invariant |
[out] | d2ibdc2_ani | 2nd derivative of the 4th anisotropic invariant |
Definition at line 143 of file Hyperelastic.f90.