DAMASK_EICMD/src/kinematics_thermal_expansio...

126 lines
5.3 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
!> @brief material subroutine incorporating kinematics resulting from thermal expansion
!> @details to be done
!--------------------------------------------------------------------------------------------------
submodule(constitutive:constitutive_thermal) kinematics_thermal_expansion
2020-02-28 14:55:07 +05:30
2020-02-29 14:06:42 +05:30
integer, dimension(:), allocatable :: kinematics_thermal_expansion_instance
2019-05-17 02:44:47 +05:30
type :: tParameters
2020-02-29 15:40:23 +05:30
real(pReal) :: &
T_ref
2020-02-29 16:51:03 +05:30
real(pReal), dimension(3,3,3) :: &
2020-09-23 04:46:12 +05:30
A = 0.0_pReal
end type tParameters
2020-02-28 14:55:07 +05:30
type(tParameters), dimension(:), allocatable :: param
2020-02-28 14:55:07 +05:30
contains
!--------------------------------------------------------------------------------------------------
!> @brief module initialization
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
2020-08-15 19:32:10 +05:30
module function kinematics_thermal_expansion_init(kinematics_length) result(myKinematics)
2020-02-28 14:55:07 +05:30
integer, intent(in) :: kinematics_length
2020-08-15 19:32:10 +05:30
logical, dimension(:,:), allocatable :: myKinematics
2020-02-28 14:55:07 +05:30
integer :: Ninstances,p,i,k
2020-08-15 19:32:10 +05:30
real(pReal), dimension(:), allocatable :: temp
class(tNode), pointer :: &
phases, &
phase, &
kinematics, &
kinematic_type
2020-09-18 02:27:56 +05:30
print'(/,a)', ' <<<+- kinematics_thermal_expansion init -+>>>'
2020-08-15 19:32:10 +05:30
myKinematics = kinematics_active('thermal_expansion',kinematics_length)
Ninstances = count(myKinematics)
print'(a,i2)', ' # instances: ',Ninstances; flush(IO_STDOUT)
if(Ninstances == 0) return
2020-02-28 14:55:07 +05:30
phases => config_material%get('phase')
allocate(param(Ninstances))
2020-08-15 19:32:10 +05:30
allocate(kinematics_thermal_expansion_instance(phases%length), source=0)
do p = 1, phases%length
if(any(myKinematics(:,p))) kinematics_thermal_expansion_instance(p) = count(myKinematics(:,1:p))
phase => phases%get(p)
2020-08-15 19:32:10 +05:30
if(count(myKinematics(:,p)) == 0) cycle
kinematics => phase%get('kinematics')
do k = 1, kinematics%length
if(myKinematics(k,p)) then
associate(prm => param(kinematics_thermal_expansion_instance(p)))
kinematic_type => kinematics%get(k)
2020-08-15 19:32:10 +05:30
prm%T_ref = kinematic_type%get_asFloat('T_ref', defaultVal=0.0_pReal)
! read up to three parameters (constant, linear, quadratic with T)
temp = kinematic_type%get_asFloats('A_11')
2020-09-23 04:46:12 +05:30
prm%A(1,1,1:size(temp)) = temp
2020-08-15 19:32:10 +05:30
temp = kinematic_type%get_asFloats('A_22',defaultVal=[(0.0_pReal, i=1,size(temp))],requiredSize=size(temp))
2020-09-23 04:46:12 +05:30
prm%A(2,2,1:size(temp)) = temp
2020-08-15 19:32:10 +05:30
temp = kinematic_type%get_asFloats('A_33',defaultVal=[(0.0_pReal, i=1,size(temp))],requiredSize=size(temp))
2020-09-23 04:46:12 +05:30
prm%A(3,3,1:size(temp)) = temp
do i=1, size(prm%A,3)
prm%A(1:3,1:3,i) = lattice_applyLatticeSymmetry33(prm%A(1:3,1:3,i),&
2020-08-15 19:32:10 +05:30
phase%get_asString('lattice'))
enddo
end associate
endif
2020-02-29 17:27:19 +05:30
enddo
enddo
2020-08-15 19:32:10 +05:30
end function kinematics_thermal_expansion_init
!--------------------------------------------------------------------------------------------------
2020-06-26 15:14:17 +05:30
!> @brief constitutive equation for calculating the velocity gradient
!--------------------------------------------------------------------------------------------------
2020-12-23 11:37:18 +05:30
module subroutine kinematics_thermal_expansion_LiAndItsTangent(Li, dLi_dTstar, co, ip, el)
2020-02-28 14:55:07 +05:30
2019-05-17 02:44:47 +05:30
integer, intent(in) :: &
2020-12-23 11:37:18 +05:30
co, & !< grain number
2019-06-16 00:02:53 +05:30
ip, & !< integration point number
el !< element number
real(pReal), intent(out), dimension(3,3) :: &
2019-06-16 00:02:53 +05:30
Li !< thermal velocity gradient
real(pReal), intent(out), dimension(3,3,3,3) :: &
2019-06-16 00:02:53 +05:30
dLi_dTstar !< derivative of Li with respect to Tstar (4th-order tensor defined to be zero)
2020-02-29 14:06:42 +05:30
2019-05-17 02:44:47 +05:30
integer :: &
phase, &
2020-02-29 14:06:42 +05:30
homog
real(pReal) :: &
2020-02-29 15:40:23 +05:30
T, TDot
2020-02-28 14:55:07 +05:30
2020-12-23 11:37:18 +05:30
phase = material_phaseAt(co,el)
2019-03-10 15:32:32 +05:30
homog = material_homogenizationAt(el)
T = temperature(homog)%p(material_homogenizationMemberAt(ip,el))
TDot = temperatureRate(homog)%p(material_homogenizationMemberAt(ip,el))
2020-02-28 14:55:07 +05:30
2020-02-29 15:40:23 +05:30
associate(prm => param(kinematics_thermal_expansion_instance(phase)))
Li = TDot * ( &
2020-09-23 04:46:12 +05:30
prm%A(1:3,1:3,1)*(T - prm%T_ref)**0 & ! constant coefficient
+ prm%A(1:3,1:3,2)*(T - prm%T_ref)**1 & ! linear coefficient
+ prm%A(1:3,1:3,3)*(T - prm%T_ref)**2 & ! quadratic coefficient
) / &
(1.0_pReal &
2020-09-23 04:46:12 +05:30
+ prm%A(1:3,1:3,1)*(T - prm%T_ref)**1 / 1. &
+ prm%A(1:3,1:3,2)*(T - prm%T_ref)**2 / 2. &
+ prm%A(1:3,1:3,3)*(T - prm%T_ref)**3 / 3. &
)
2020-02-29 15:40:23 +05:30
end associate
2020-02-28 14:55:07 +05:30
dLi_dTstar = 0.0_pReal
end subroutine kinematics_thermal_expansion_LiAndItsTangent
end submodule kinematics_thermal_expansion