2015-05-28 22:32:23 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
|
|
|
|
!> @brief material subroutine incorporating kinematics resulting from thermal expansion
|
|
|
|
!> @details to be done
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
module kinematics_thermal_expansion
|
2019-05-15 02:42:32 +05:30
|
|
|
use prec
|
2019-05-17 02:44:47 +05:30
|
|
|
use IO
|
|
|
|
use config
|
|
|
|
use debug
|
|
|
|
use math
|
|
|
|
use lattice
|
|
|
|
use material
|
|
|
|
|
2019-02-23 01:07:41 +05:30
|
|
|
implicit none
|
|
|
|
private
|
2019-02-13 04:05:22 +05:30
|
|
|
|
2019-05-17 02:44:47 +05:30
|
|
|
type :: tParameters
|
2019-02-23 01:07:41 +05:30
|
|
|
real(pReal), allocatable, dimension(:,:,:) :: &
|
|
|
|
expansion
|
|
|
|
end type tParameters
|
2019-02-13 04:05:22 +05:30
|
|
|
|
2019-02-23 01:07:41 +05:30
|
|
|
type(tParameters), dimension(:), allocatable :: param
|
|
|
|
|
|
|
|
public :: &
|
|
|
|
kinematics_thermal_expansion_init, &
|
|
|
|
kinematics_thermal_expansion_initialStrain, &
|
|
|
|
kinematics_thermal_expansion_LiAndItsTangent
|
2015-05-28 22:32:23 +05:30
|
|
|
|
|
|
|
contains
|
|
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief module initialization
|
|
|
|
!> @details reads in material parameters, allocates arrays, and does sanity checks
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-05-17 02:44:47 +05:30
|
|
|
subroutine kinematics_thermal_expansion_init
|
2018-12-31 01:28:38 +05:30
|
|
|
|
2019-05-17 02:44:47 +05:30
|
|
|
integer :: &
|
2019-02-23 01:07:41 +05:30
|
|
|
Ninstance, &
|
|
|
|
p, i
|
|
|
|
real(pReal), dimension(:), allocatable :: &
|
|
|
|
temp
|
|
|
|
|
|
|
|
write(6,'(/,a)') ' <<<+- kinematics_'//KINEMATICS_thermal_expansion_LABEL//' init -+>>>'
|
2015-05-28 22:32:23 +05:30
|
|
|
|
2019-05-17 02:44:47 +05:30
|
|
|
Ninstance = count(phase_kinematics == KINEMATICS_thermal_expansion_ID)
|
2019-02-23 01:07:41 +05:30
|
|
|
|
2019-05-17 02:44:47 +05:30
|
|
|
if (iand(debug_level(debug_constitutive),debug_levelBasic) /= 0) &
|
2019-02-23 01:07:41 +05:30
|
|
|
write(6,'(a16,1x,i5,/)') '# instances:',Ninstance
|
|
|
|
|
|
|
|
allocate(param(Ninstance))
|
|
|
|
|
2019-05-17 02:44:47 +05:30
|
|
|
do p = 1, size(phase_kinematics)
|
2019-02-23 01:07:41 +05:30
|
|
|
if (all(phase_kinematics(:,p) /= KINEMATICS_thermal_expansion_ID)) cycle
|
|
|
|
|
|
|
|
! ToDo: Here we need to decide how to extend the concept of instances to
|
|
|
|
! kinetics and sources. I would suggest that the same mechanism exists at maximum once per phase
|
|
|
|
|
|
|
|
! read up to three parameters (constant, linear, quadratic with T)
|
|
|
|
temp = config_phase(p)%getFloats('thermal_expansion11')
|
|
|
|
!lattice_thermalExpansion33(1,1,1:size(temp),p) = temp
|
|
|
|
temp = config_phase(p)%getFloats('thermal_expansion22', &
|
|
|
|
defaultVal=[(0.0_pReal, i=1,size(temp))],requiredSize=size(temp))
|
|
|
|
!lattice_thermalExpansion33(2,2,1:size(temp),p) = temp
|
|
|
|
temp = config_phase(p)%getFloats('thermal_expansion33', &
|
|
|
|
defaultVal=[(0.0_pReal, i=1,size(temp))],requiredSize=size(temp))
|
|
|
|
enddo
|
2015-05-28 22:32:23 +05:30
|
|
|
|
|
|
|
end subroutine kinematics_thermal_expansion_init
|
|
|
|
|
2019-02-23 01:07:41 +05:30
|
|
|
|
2015-07-23 21:29:25 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief report initial thermal strain based on current temperature deviation from reference
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-02-13 04:05:22 +05:30
|
|
|
pure function kinematics_thermal_expansion_initialStrain(homog,phase,offset)
|
2019-02-23 01:07:41 +05:30
|
|
|
|
2019-05-17 02:44:47 +05:30
|
|
|
integer, intent(in) :: &
|
2019-02-23 01:07:41 +05:30
|
|
|
phase, &
|
|
|
|
homog, offset
|
|
|
|
real(pReal), dimension(3,3) :: &
|
|
|
|
kinematics_thermal_expansion_initialStrain !< initial thermal strain (should be small strain, though)
|
2015-07-23 21:29:25 +05:30
|
|
|
|
2019-02-23 01:07:41 +05:30
|
|
|
|
|
|
|
kinematics_thermal_expansion_initialStrain = &
|
|
|
|
(temperature(homog)%p(offset) - lattice_referenceTemperature(phase))**1 / 1. * &
|
|
|
|
lattice_thermalExpansion33(1:3,1:3,1,phase) + & ! constant coefficient
|
|
|
|
(temperature(homog)%p(offset) - lattice_referenceTemperature(phase))**2 / 2. * &
|
|
|
|
lattice_thermalExpansion33(1:3,1:3,2,phase) + & ! linear coefficient
|
|
|
|
(temperature(homog)%p(offset) - lattice_referenceTemperature(phase))**3 / 3. * &
|
|
|
|
lattice_thermalExpansion33(1:3,1:3,3,phase) ! quadratic coefficient
|
2015-07-23 21:29:25 +05:30
|
|
|
|
2015-07-24 20:17:18 +05:30
|
|
|
end function kinematics_thermal_expansion_initialStrain
|
2015-07-23 21:29:25 +05:30
|
|
|
|
2019-02-23 01:07:41 +05:30
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief contains the constitutive equation for calculating the velocity gradient
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-12-31 01:28:38 +05:30
|
|
|
subroutine kinematics_thermal_expansion_LiAndItsTangent(Li, dLi_dTstar, ipc, ip, el)
|
2019-02-23 01:07:41 +05:30
|
|
|
|
2019-05-17 02:44:47 +05:30
|
|
|
integer, intent(in) :: &
|
2019-02-23 01:07:41 +05:30
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
|
|
|
real(pReal), intent(out), dimension(3,3) :: &
|
|
|
|
Li !< thermal velocity gradient
|
|
|
|
real(pReal), intent(out), dimension(3,3,3,3) :: &
|
|
|
|
dLi_dTstar !< derivative of Li with respect to Tstar (4th-order tensor defined to be zero)
|
2019-05-17 02:44:47 +05:30
|
|
|
integer :: &
|
2019-02-23 01:07:41 +05:30
|
|
|
phase, &
|
|
|
|
homog, offset
|
|
|
|
real(pReal) :: &
|
|
|
|
T, TRef, TDot
|
|
|
|
|
|
|
|
phase = material_phase(ipc,ip,el)
|
2019-03-10 15:32:32 +05:30
|
|
|
homog = material_homogenizationAt(el)
|
2019-02-23 01:07:41 +05:30
|
|
|
offset = thermalMapping(homog)%p(ip,el)
|
|
|
|
T = temperature(homog)%p(offset)
|
|
|
|
TDot = temperatureRate(homog)%p(offset)
|
|
|
|
TRef = lattice_referenceTemperature(phase)
|
|
|
|
|
|
|
|
Li = TDot * ( &
|
|
|
|
lattice_thermalExpansion33(1:3,1:3,1,phase)*(T - TRef)**0 & ! constant coefficient
|
|
|
|
+ lattice_thermalExpansion33(1:3,1:3,2,phase)*(T - TRef)**1 & ! linear coefficient
|
|
|
|
+ lattice_thermalExpansion33(1:3,1:3,3,phase)*(T - TRef)**2 & ! quadratic coefficient
|
|
|
|
) / &
|
|
|
|
(1.0_pReal &
|
|
|
|
+ lattice_thermalExpansion33(1:3,1:3,1,phase)*(T - TRef)**1 / 1. &
|
|
|
|
+ lattice_thermalExpansion33(1:3,1:3,2,phase)*(T - TRef)**2 / 2. &
|
|
|
|
+ lattice_thermalExpansion33(1:3,1:3,3,phase)*(T - TRef)**3 / 3. &
|
|
|
|
)
|
|
|
|
dLi_dTstar = 0.0_pReal
|
2015-05-28 22:32:23 +05:30
|
|
|
|
|
|
|
end subroutine kinematics_thermal_expansion_LiAndItsTangent
|
|
|
|
|
|
|
|
end module kinematics_thermal_expansion
|