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
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-07-09 04:31:08 +05:30
|
|
|
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
|
2019-02-23 01:07:41 +05:30
|
|
|
end type tParameters
|
2020-02-28 14:55:07 +05:30
|
|
|
|
2019-02-23 01:07:41 +05:30
|
|
|
type(tParameters), dimension(:), allocatable :: param
|
2020-02-28 14:55:07 +05:30
|
|
|
|
2015-05-28 22:32:23 +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
|
|
|
|
2020-08-15 19:32:10 +05:30
|
|
|
integer, intent(in) :: kinematics_length
|
|
|
|
logical, dimension(:,:), allocatable :: myKinematics
|
2020-02-28 14:55:07 +05:30
|
|
|
|
2020-08-15 19:32:10 +05:30
|
|
|
integer :: Ninstance,p,i,k
|
|
|
|
real(pReal), dimension(:), allocatable :: temp
|
|
|
|
class(tNode), pointer :: &
|
|
|
|
phases, &
|
|
|
|
phase, &
|
|
|
|
pl, &
|
|
|
|
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)
|
|
|
|
Ninstance = count(myKinematics)
|
2020-09-22 16:39:12 +05:30
|
|
|
print'(a,i2)', ' # instances: ',Ninstance; flush(IO_STDOUT)
|
2020-08-15 19:32:10 +05:30
|
|
|
if(Ninstance == 0) return
|
2020-02-28 14:55:07 +05:30
|
|
|
|
2020-09-13 14:09:17 +05:30
|
|
|
phases => config_material%get('phase')
|
2019-02-23 01:07:41 +05:30
|
|
|
allocate(param(Ninstance))
|
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)
|
|
|
|
pl => phase%get('plasticity')
|
|
|
|
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)
|
|
|
|
|
|
|
|
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
|
2019-02-23 01:07:41 +05:30
|
|
|
enddo
|
2015-05-28 22:32:23 +05:30
|
|
|
|
2020-08-15 19:32:10 +05:30
|
|
|
|
|
|
|
end function kinematics_thermal_expansion_init
|
2015-05-28 22:32:23 +05:30
|
|
|
|
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
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-07-09 05:19:48 +05:30
|
|
|
pure module function kinematics_thermal_expansion_initialStrain(homog,phase,offset) result(initialStrain)
|
|
|
|
|
|
|
|
integer, intent(in) :: &
|
|
|
|
phase, &
|
|
|
|
homog, &
|
|
|
|
offset
|
|
|
|
|
|
|
|
real(pReal), dimension(3,3) :: &
|
2020-09-23 04:46:12 +05:30
|
|
|
initialStrain !< initial thermal strain (should be small strain, though)
|
2020-07-09 05:19:48 +05:30
|
|
|
|
|
|
|
associate(prm => param(kinematics_thermal_expansion_instance(phase)))
|
|
|
|
initialStrain = &
|
2020-09-23 04:46:12 +05:30
|
|
|
(temperature(homog)%p(offset) - prm%T_ref)**1 / 1. * prm%A(1:3,1:3,1) + & ! constant coefficient
|
|
|
|
(temperature(homog)%p(offset) - prm%T_ref)**2 / 2. * prm%A(1:3,1:3,2) + & ! linear coefficient
|
|
|
|
(temperature(homog)%p(offset) - prm%T_ref)**3 / 3. * prm%A(1:3,1:3,3) ! quadratic coefficient
|
2020-07-09 05:19:48 +05:30
|
|
|
end associate
|
|
|
|
|
|
|
|
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
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-06-26 15:14:17 +05:30
|
|
|
!> @brief constitutive equation for calculating the velocity gradient
|
2015-05-28 22:32:23 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-07-09 04:31:08 +05:30
|
|
|
module subroutine kinematics_thermal_expansion_LiAndItsTangent(Li, dLi_dTstar, ipc, ip, el)
|
2020-02-28 14:55:07 +05:30
|
|
|
|
2019-05-17 02:44:47 +05:30
|
|
|
integer, intent(in) :: &
|
2019-06-16 00:02:53 +05:30
|
|
|
ipc, & !< grain number
|
|
|
|
ip, & !< integration point number
|
|
|
|
el !< element number
|
2019-02-23 01:07:41 +05:30
|
|
|
real(pReal), intent(out), dimension(3,3) :: &
|
2019-06-16 00:02:53 +05:30
|
|
|
Li !< thermal velocity gradient
|
2019-02-23 01:07:41 +05:30
|
|
|
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 :: &
|
2019-02-23 01:07:41 +05:30
|
|
|
phase, &
|
2020-02-29 14:06:42 +05:30
|
|
|
homog
|
2019-02-23 01:07:41 +05:30
|
|
|
real(pReal) :: &
|
2020-02-29 15:40:23 +05:30
|
|
|
T, TDot
|
2020-02-28 14:55:07 +05:30
|
|
|
|
2019-06-14 13:07:01 +05:30
|
|
|
phase = material_phaseAt(ipc,el)
|
2019-03-10 15:32:32 +05:30
|
|
|
homog = material_homogenizationAt(el)
|
2020-02-29 14:06:42 +05:30
|
|
|
T = temperature(homog)%p(thermalMapping(homog)%p(ip,el))
|
|
|
|
TDot = temperatureRate(homog)%p(thermalMapping(homog)%p(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)))
|
2019-02-23 01:07:41 +05:30
|
|
|
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
|
2019-02-23 01:07:41 +05:30
|
|
|
) / &
|
|
|
|
(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. &
|
2019-02-23 01:07:41 +05:30
|
|
|
)
|
2020-02-29 15:40:23 +05:30
|
|
|
end associate
|
2020-02-28 14:55:07 +05:30
|
|
|
dLi_dTstar = 0.0_pReal
|
|
|
|
|
2015-05-28 22:32:23 +05:30
|
|
|
end subroutine kinematics_thermal_expansion_LiAndItsTangent
|
|
|
|
|
2020-07-09 04:31:08 +05:30
|
|
|
end submodule kinematics_thermal_expansion
|