This commit is contained in:
Philip Eisenlohr 2021-11-21 15:49:04 -05:00
parent 021d614daf
commit da23c916ca
4 changed files with 105 additions and 107 deletions

View File

@ -405,7 +405,7 @@ module lattice
contains contains
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Module initialization !> @brief module initialization
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
subroutine lattice_init subroutine lattice_init
@ -417,7 +417,7 @@ end subroutine lattice_init
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Characteristic shear for twinning !> @brief characteristic shear for twinning
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function lattice_characteristicShear_Twin(Ntwin,lattice,CoverA) result(characteristicShear) function lattice_characteristicShear_Twin(Ntwin,lattice,CoverA) result(characteristicShear)
@ -491,7 +491,7 @@ end function lattice_characteristicShear_Twin
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Rotated elasticity matrices for twinning in 66-vector notation !> @brief rotated elasticity matrices for twinning in 6x6-matrix notation
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function lattice_C66_twin(Ntwin,C66,lattice,CoverA) function lattice_C66_twin(Ntwin,C66,lattice,CoverA)
@ -529,7 +529,7 @@ end function lattice_C66_twin
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Rotated elasticity matrices for transformation in 66-vector notation !> @brief rotated elasticity matrices for transformation in 6x6-matrix notation
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function lattice_C66_trans(Ntrans,C_parent66,lattice_target, & function lattice_C66_trans(Ntrans,C_parent66,lattice_target, &
cOverA_trans,a_bcc,a_fcc) cOverA_trans,a_bcc,a_fcc)
@ -588,7 +588,7 @@ function lattice_C66_trans(Ntrans,C_parent66,lattice_target, &
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Non-schmid projections for bcc with up to 6 coefficients !> @brief non-Schmid projections for bcc with up to 6 coefficients
! Koester et al. 2012, Acta Materialia 60 (2012) 38943901, eq. (17) ! Koester et al. 2012, Acta Materialia 60 (2012) 38943901, eq. (17)
! Gröger et al. 2008, Acta Materialia 56 (2008) 54125425, table 1 ! Gröger et al. 2008, Acta Materialia 56 (2008) 54125425, table 1
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
@ -635,7 +635,7 @@ end function lattice_nonSchmidMatrix
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Slip-slip interaction matrix !> @brief slip-slip interaction matrix
!> details only active slip systems are considered !> details only active slip systems are considered
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function lattice_interaction_SlipBySlip(Nslip,interactionValues,lattice) result(interactionMatrix) function lattice_interaction_SlipBySlip(Nslip,interactionValues,lattice) result(interactionMatrix)
@ -883,7 +883,7 @@ end function lattice_interaction_SlipBySlip
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Twin-twin interaction matrix !> @brief twin-twin interaction matrix
!> details only active twin systems are considered !> details only active twin systems are considered
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function lattice_interaction_TwinByTwin(Ntwin,interactionValues,lattice) result(interactionMatrix) function lattice_interaction_TwinByTwin(Ntwin,interactionValues,lattice) result(interactionMatrix)
@ -981,7 +981,7 @@ end function lattice_interaction_TwinByTwin
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Trans-trans interaction matrix !> @brief trans-trans interaction matrix
!> details only active trans systems are considered !> details only active trans systems are considered
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function lattice_interaction_TransByTrans(Ntrans,interactionValues,lattice) result(interactionMatrix) function lattice_interaction_TransByTrans(Ntrans,interactionValues,lattice) result(interactionMatrix)
@ -1023,7 +1023,7 @@ end function lattice_interaction_TransByTrans
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Slip-twin interaction matrix !> @brief slip-twin interaction matrix
!> details only active slip and twin systems are considered !> details only active slip and twin systems are considered
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function lattice_interaction_SlipByTwin(Nslip,Ntwin,interactionValues,lattice) result(interactionMatrix) function lattice_interaction_SlipByTwin(Nslip,Ntwin,interactionValues,lattice) result(interactionMatrix)
@ -1186,7 +1186,7 @@ end function lattice_interaction_SlipByTwin
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Slip-trans interaction matrix !> @brief slip-trans interaction matrix
!> details only active slip and trans systems are considered !> details only active slip and trans systems are considered
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function lattice_interaction_SlipByTrans(Nslip,Ntrans,interactionValues,lattice) result(interactionMatrix) function lattice_interaction_SlipByTrans(Nslip,Ntrans,interactionValues,lattice) result(interactionMatrix)
@ -1239,7 +1239,7 @@ function lattice_interaction_SlipByTrans(Nslip,Ntrans,interactionValues,lattice)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Twin-slip interaction matrix !> @brief twin-slip interaction matrix
!> details only active twin and slip systems are considered !> details only active twin and slip systems are considered
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function lattice_interaction_TwinBySlip(Ntwin,Nslip,interactionValues,lattice) result(interactionMatrix) function lattice_interaction_TwinBySlip(Ntwin,Nslip,interactionValues,lattice) result(interactionMatrix)
@ -1411,7 +1411,7 @@ end function lattice_SchmidMatrix_twin
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Schmid matrix for twinning !> @brief Schmid matrix for transformation
!> details only active twin systems are considered !> details only active twin systems are considered
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function lattice_SchmidMatrix_trans(Ntrans,lattice_target,cOverA,a_bcc,a_fcc) result(SchmidMatrix) function lattice_SchmidMatrix_trans(Ntrans,lattice_target,cOverA,a_bcc,a_fcc) result(SchmidMatrix)
@ -1483,7 +1483,7 @@ end function lattice_SchmidMatrix_cleavage
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Slip direction of slip systems (|| b) !> @brief slip direction of slip systems (|| b)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function lattice_slip_direction(Nslip,lattice,cOverA) result(d) function lattice_slip_direction(Nslip,lattice,cOverA) result(d)
@ -1501,7 +1501,7 @@ end function lattice_slip_direction
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Normal direction of slip systems (|| n) !> @brief normal direction of slip systems (|| n)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function lattice_slip_normal(Nslip,lattice,cOverA) result(n) function lattice_slip_normal(Nslip,lattice,cOverA) result(n)
@ -1519,7 +1519,7 @@ end function lattice_slip_normal
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Transverse direction of slip systems (|| t = b x n) !> @brief transverse direction of slip systems (|| t = b x n)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function lattice_slip_transverse(Nslip,lattice,cOverA) result(t) function lattice_slip_transverse(Nslip,lattice,cOverA) result(t)
@ -1537,7 +1537,7 @@ end function lattice_slip_transverse
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Labels for slip systems !> @brief labels of slip systems
!> details only active slip systems are considered !> details only active slip systems are considered
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function lattice_labels_slip(Nslip,lattice) result(labels) function lattice_labels_slip(Nslip,lattice) result(labels)
@ -1578,7 +1578,7 @@ end function lattice_labels_slip
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Return 3x3 tensor with symmetry according to given Bravais lattice !> @brief return 3x3 tensor with symmetry according to given Bravais lattice
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
pure function lattice_symmetrize_33(T,lattice) result(T_sym) pure function lattice_symmetrize_33(T,lattice) result(T_sym)
@ -1605,7 +1605,7 @@ end function lattice_symmetrize_33
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Return stiffness matrix in 6x6 notation with symmetry according to given Bravais lattice !> @brief return stiffness matrix in 6x6 notation with symmetry according to given Bravais lattice
!> @details J. A. Rayne and B. S. Chandrasekhar Phys. Rev. 120, 1658 Erratum Phys. Rev. 122, 1962 !> @details J. A. Rayne and B. S. Chandrasekhar Phys. Rev. 120, 1658 Erratum Phys. Rev. 122, 1962
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
pure function lattice_symmetrize_C66(C66,lattice) result(C66_sym) pure function lattice_symmetrize_C66(C66,lattice) result(C66_sym)
@ -1651,7 +1651,7 @@ end function lattice_symmetrize_C66
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Labels for twin systems !> @brief labels of twin systems
!> details only active twin systems are considered !> details only active twin systems are considered
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function lattice_labels_twin(Ntwin,lattice) result(labels) function lattice_labels_twin(Ntwin,lattice) result(labels)
@ -1689,7 +1689,7 @@ end function lattice_labels_twin
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Projection of the transverse direction onto the slip plane !> @brief projection of the transverse direction onto the slip plane
!> @details: This projection is used to calculate forest hardening for edge dislocations !> @details: This projection is used to calculate forest hardening for edge dislocations
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function slipProjection_transverse(Nslip,lattice,cOverA) result(projection) function slipProjection_transverse(Nslip,lattice,cOverA) result(projection)
@ -1713,7 +1713,7 @@ end function slipProjection_transverse
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Projection of the slip direction onto the slip plane !> @brief projection of the slip direction onto the slip plane
!> @details: This projection is used to calculate forest hardening for screw dislocations !> @details: This projection is used to calculate forest hardening for screw dislocations
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function slipProjection_direction(Nslip,lattice,cOverA) result(projection) function slipProjection_direction(Nslip,lattice,cOverA) result(projection)
@ -1779,7 +1779,7 @@ end function coordinateSystem_slip
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Populate reduced interaction matrix !> @brief populate reduced interaction matrix
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function buildInteraction(reacting_used,acting_used,reacting_max,acting_max,values,matrix) function buildInteraction(reacting_used,acting_used,reacting_max,acting_max,values,matrix)
@ -1822,7 +1822,7 @@ end function buildInteraction
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Build a local coordinate system on slip, twin, trans, cleavage systems !> @brief build a local coordinate system on slip, twin, trans, cleavage systems
!> @details Order: Direction, plane (normal), and common perpendicular !> @details Order: Direction, plane (normal), and common perpendicular
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function buildCoordinateSystem(active,potential,system,lattice,cOverA) function buildCoordinateSystem(active,potential,system,lattice,cOverA)
@ -1889,7 +1889,7 @@ end function buildCoordinateSystem
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Helper function to define transformation systems !> @brief helper function to define transformation systems
! Needed to calculate Schmid matrix and rotated stiffness matrices. ! Needed to calculate Schmid matrix and rotated stiffness matrices.
! @details: set c/a = 0.0 for fcc -> bcc transformation ! @details: set c/a = 0.0 for fcc -> bcc transformation
! set a_Xcc = 0.0 for fcc -> hex transformation ! set a_Xcc = 0.0 for fcc -> hex transformation
@ -2073,7 +2073,7 @@ end function getlabels
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Equivalent Poisson's ratio (ν) !> @brief equivalent Poisson's ratio (ν)
!> @details https://doi.org/10.1143/JPSJ.20.635 !> @details https://doi.org/10.1143/JPSJ.20.635
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function lattice_equivalent_nu(C,assumption) result(nu) function lattice_equivalent_nu(C,assumption) result(nu)
@ -2106,7 +2106,7 @@ end function lattice_equivalent_nu
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Equivalent shear modulus (μ) !> @brief equivalent shear modulus (μ)
!> @details https://doi.org/10.1143/JPSJ.20.635 !> @details https://doi.org/10.1143/JPSJ.20.635
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function lattice_equivalent_mu(C,assumption) result(mu) function lattice_equivalent_mu(C,assumption) result(mu)
@ -2135,7 +2135,7 @@ end function lattice_equivalent_mu
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Check correctness of some lattice functions. !> @brief check correctness of some lattice functions
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
subroutine selfTest subroutine selfTest

View File

@ -15,15 +15,15 @@ module math
implicit none implicit none
public public
#if __INTEL_COMPILER >= 1900 #if __INTEL_COMPILER >= 1900
! do not make use associated entities available to other modules ! do not make use of associated entities available to other modules
private :: & private :: &
IO, & IO, &
config config
#endif #endif
real(pReal), parameter :: PI = acos(-1.0_pReal) !< ratio of a circle's circumference to its diameter real(pReal), parameter :: PI = acos(-1.0_pReal) !< ratio of a circle's circumference to its diameter
real(pReal), parameter :: INDEG = 180.0_pReal/PI !< conversion from radian into degree real(pReal), parameter :: INDEG = 180.0_pReal/PI !< conversion from radian to degree
real(pReal), parameter :: INRAD = PI/180.0_pReal !< conversion from degree into radian real(pReal), parameter :: INRAD = PI/180.0_pReal !< conversion from degree to radian
complex(pReal), parameter :: TWOPIIMG = cmplx(0.0_pReal,2.0_pReal*PI) !< Re(0.0), Im(2xPi) complex(pReal), parameter :: TWOPIIMG = cmplx(0.0_pReal,2.0_pReal*PI) !< Re(0.0), Im(2xPi)
real(pReal), dimension(3,3), parameter :: & real(pReal), dimension(3,3), parameter :: &
@ -822,7 +822,7 @@ end function math_sym3333to66
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief convert 66 matrix into symmetric 3x3x3x3 matrix !> @brief convert 6x6 matrix into symmetric 3x3x3x3 matrix
!> @details Weighted conversion (default) rearranges according to Nye and weights shear !> @details Weighted conversion (default) rearranges according to Nye and weights shear
! components according to Mandel. Advisable for matrix operations. ! components according to Mandel. Advisable for matrix operations.
! Unweighted conversion only rearranges order according to Nye ! Unweighted conversion only rearranges order according to Nye
@ -854,7 +854,7 @@ end function math_66toSym3333
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Convert 6x6 Voigt matrix into symmetric 3x3x3x3 matrix. !> @brief convert 6x6 Voigt matrix into symmetric 3x3x3x3 matrix
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
pure function math_Voigt66to3333(m66) pure function math_Voigt66to3333(m66)
@ -875,7 +875,7 @@ end function math_Voigt66to3333
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Convert symmetric 3x3x3x3 matrix into 6x6 Voigt matrix. !> @brief convert symmetric 3x3x3x3 matrix into 6x6 Voigt matrix
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
pure function math_3333toVoigt66(m3333) pure function math_3333toVoigt66(m3333)

View File

@ -15,7 +15,7 @@ submodule(phase:mechanical) elastic
contains contains
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Initialize elasticity. !> @brief initialize elasticity
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
module subroutine elastic_init(phases) module subroutine elastic_init(phases)
@ -62,7 +62,7 @@ end subroutine elastic_init
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Return 6x6 elasticity tensor. !> @brief return 6x6 elasticity tensor
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
module function elastic_C66(ph,en) result(C66) module function elastic_C66(ph,en) result(C66)
@ -93,7 +93,7 @@ end function elastic_C66
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Return shear modulus. !> @brief return shear modulus
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
module function elastic_mu(ph,en) result(mu) module function elastic_mu(ph,en) result(mu)
@ -110,7 +110,7 @@ end function elastic_mu
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Return Poisson ratio. !> @brief return Poisson ratio
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
module function elastic_nu(ph,en) result(nu) module function elastic_nu(ph,en) result(nu)
@ -128,7 +128,7 @@ end function elastic_nu
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief returns the 2nd Piola-Kirchhoff stress tensor and its tangent with respect to !> @brief return the 2nd Piola-Kirchhoff stress tensor and its tangent with respect to
!> the elastic and intermediate deformation gradients using Hooke's law !> the elastic and intermediate deformation gradients using Hooke's law
! ToDo: Use Voigt matrix directly ! ToDo: Use Voigt matrix directly
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------

View File

@ -473,9 +473,7 @@ module function plastic_dislotwin_homogenizedC(ph,en) result(homogenizedC)
integer, intent(in) :: & integer, intent(in) :: &
ph, en ph, en
real(pReal), dimension(6,6) :: & real(pReal), dimension(6,6) :: &
homogenizedC homogenizedC, &
real(pReal), dimension(6,6) :: &
C C
real(pReal), dimension(:,:,:), allocatable :: & real(pReal), dimension(:,:,:), allocatable :: &
C66_tw, & C66_tw, &