DAMASK_EICMD/src/homogenization.f90

483 lines
17 KiB
Fortran
Raw Normal View History

2012-08-10 21:28:17 +05:30
!--------------------------------------------------------------------------------------------------
!> @author Franz Roters, Max-Planck-Institut für Eisenforschung GmbH
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
!> @author Denny Tjahjanto, Max-Planck-Institut für Eisenforschung GmbH
!> @brief homogenization manager, organizing deformation partitioning and stress homogenization
2012-08-10 21:28:17 +05:30
!--------------------------------------------------------------------------------------------------
module homogenization
2019-06-15 20:12:16 +05:30
use prec
use IO
use config
use math
use material
use phase
2019-06-15 20:12:16 +05:30
use discretization
2020-12-30 04:44:48 +05:30
use HDF5_utilities
2019-06-15 20:12:16 +05:30
use results
use lattice
2020-03-17 05:09:32 +05:30
2019-06-15 20:12:16 +05:30
implicit none
private
2019-06-11 13:18:07 +05:30
2021-02-12 20:01:43 +05:30
enum, bind(c); enumerator :: &
THERMAL_ISOTHERMAL_ID, &
THERMAL_CONDUCTION_ID, &
DAMAGE_NONE_ID, &
DAMAGE_NONLOCAL_ID, &
HOMOGENIZATION_UNDEFINED_ID, &
HOMOGENIZATION_NONE_ID, &
HOMOGENIZATION_ISOSTRAIN_ID, &
HOMOGENIZATION_RGC_ID
end enum
type(tState), allocatable, dimension(:), public :: &
homogState, &
damageState_h
integer(kind(THERMAL_isothermal_ID)), dimension(:), allocatable :: &
2021-02-12 20:01:43 +05:30
thermal_type !< thermal transport model
2021-03-01 10:46:16 +05:30
integer(kind(DAMAGE_none_ID)), dimension(:), allocatable :: &
2021-02-12 20:01:43 +05:30
damage_type !< nonlocal damage model
type, private :: tNumerics_damage
real(pReal) :: &
charLength !< characteristic length scale for gradient problems
end type tNumerics_damage
type(tNumerics_damage), private :: &
num_damage
2020-03-29 23:34:51 +05:30
logical, public :: &
terminallyIll = .false. !< at least one material point is terminally ill
2020-04-14 13:08:32 +05:30
!--------------------------------------------------------------------------------------------------
! General variables for the homogenization at a material point
real(pReal), dimension(:,:,:), allocatable, public :: &
2020-12-16 15:51:24 +05:30
homogenization_F0, & !< def grad of IP at start of FE increment
homogenization_F !< def grad of IP to be reached at end of FE increment
real(pReal), dimension(:,:,:), allocatable, public :: & !, protected :: & Issue with ifort
2020-12-16 15:51:24 +05:30
homogenization_P !< first P--K stress of IP
real(pReal), dimension(:,:,:,:,:), allocatable, public :: & !, protected :: &
2020-12-16 15:51:24 +05:30
homogenization_dPdF !< tangent of first P--K stress at IP
2019-06-15 20:12:16 +05:30
2020-12-16 15:51:24 +05:30
!--------------------------------------------------------------------------------------------------
2020-03-17 05:09:32 +05:30
type :: tNumerics
integer :: &
2020-04-14 13:08:32 +05:30
nMPstate !< materialpoint state loop limit
2020-03-17 05:09:32 +05:30
end type tNumerics
type(tNumerics) :: num
2020-12-16 15:51:24 +05:30
!--------------------------------------------------------------------------------------------------
interface
2020-03-17 05:09:32 +05:30
2021-02-09 03:51:53 +05:30
module subroutine mechanical_init(num_homog)
2020-06-16 22:17:19 +05:30
class(tNode), pointer, intent(in) :: &
2020-12-16 15:51:24 +05:30
num_homog !< pointer to mechanical homogenization numerics data
2021-02-09 03:51:53 +05:30
end subroutine mechanical_init
2020-12-16 15:51:24 +05:30
2020-12-30 16:30:47 +05:30
module subroutine thermal_init
end subroutine thermal_init
2021-01-21 01:24:31 +05:30
module subroutine damage_init
end subroutine damage_init
2021-02-15 23:13:51 +05:30
module subroutine mechanical_partition(subF,ce)
2020-12-16 15:51:24 +05:30
real(pReal), intent(in), dimension(3,3) :: &
subF
integer, intent(in) :: &
ce
2021-02-09 03:51:53 +05:30
end subroutine mechanical_partition
2020-03-17 05:09:32 +05:30
2021-01-24 17:56:01 +05:30
module subroutine thermal_partition(ce)
integer, intent(in) :: ce
2021-01-08 02:45:18 +05:30
end subroutine thermal_partition
module subroutine damage_partition(ce)
2021-01-21 01:24:31 +05:30
integer, intent(in) :: ce
end subroutine damage_partition
2021-02-26 02:12:40 +05:30
module subroutine mechanical_homogenize(dt,ce)
real(pReal), intent(in) :: dt
2020-12-16 15:51:24 +05:30
integer, intent(in) :: &
2021-02-26 02:12:40 +05:30
ce !< cell
2021-02-09 03:51:53 +05:30
end subroutine mechanical_homogenize
2020-03-17 05:09:32 +05:30
module subroutine mechanical_results(group_base,ho)
2020-12-16 15:51:24 +05:30
character(len=*), intent(in) :: group_base
integer, intent(in) :: ho
2021-02-09 03:51:53 +05:30
end subroutine mechanical_results
2020-03-17 05:09:32 +05:30
2021-04-07 10:56:54 +05:30
module subroutine damage_results(ho,group)
integer, intent(in) :: ho
character(len=*), intent(in) :: group
end subroutine damage_results
module subroutine thermal_results(ho,group)
integer, intent(in) :: ho
character(len=*), intent(in) :: group
end subroutine thermal_results
module function mechanical_updateState(subdt,subF,ce) result(doneAndHappy)
2020-12-28 14:25:54 +05:30
real(pReal), intent(in) :: &
2021-02-22 20:47:32 +05:30
subdt !< current time step
2020-12-28 14:25:54 +05:30
real(pReal), intent(in), dimension(3,3) :: &
subF
integer, intent(in) :: &
ce !< cell
2020-12-28 14:25:54 +05:30
logical, dimension(2) :: doneAndHappy
2021-02-09 03:51:53 +05:30
end function mechanical_updateState
2019-06-15 20:12:16 +05:30
2021-04-11 12:02:13 +05:30
module function homogenization_mu_T(ce) result(mu)
integer, intent(in) :: ce
real(pReal) :: mu
end function homogenization_mu_T
2021-01-24 17:56:01 +05:30
2021-04-09 03:10:20 +05:30
module function homogenization_K_T(ce) result(K)
integer, intent(in) :: ce
2021-01-24 17:56:01 +05:30
real(pReal), dimension(3,3) :: K
2021-04-09 03:10:20 +05:30
end function homogenization_K_T
2021-01-24 17:56:01 +05:30
2021-04-11 12:02:13 +05:30
module function homogenization_f_T(ce) result(f)
2021-01-24 17:56:01 +05:30
integer, intent(in) :: ce
2021-04-11 12:02:13 +05:30
real(pReal) :: f
end function homogenization_f_T
2021-01-24 17:56:01 +05:30
module subroutine homogenization_thermal_setField(T,dot_T, ce)
integer, intent(in) :: ce
real(pReal), intent(in) :: T, dot_T
end subroutine homogenization_thermal_setField
2021-04-09 03:10:20 +05:30
module function homogenization_mu_phi(ce) result(mu)
integer, intent(in) :: ce
2021-04-09 03:10:20 +05:30
real(pReal) :: mu
2021-04-08 17:01:21 +05:30
end function homogenization_mu_phi
2021-01-24 23:17:19 +05:30
2021-04-11 11:18:54 +05:30
module function homogenization_K_phi(ce) result(K)
integer, intent(in) :: ce
real(pReal), dimension(3,3) :: K
end function homogenization_K_phi
2021-04-09 03:10:20 +05:30
module function homogenization_f_phi(phi,ce) result(f)
integer, intent(in) :: ce
2021-04-08 17:01:21 +05:30
real(pReal), intent(in) :: phi
2021-04-09 03:10:20 +05:30
real(pReal) :: f
2021-04-08 17:01:21 +05:30
end function homogenization_f_phi
2021-01-25 19:43:17 +05:30
2021-04-08 00:36:29 +05:30
module subroutine homogenization_set_phi(phi,ce)
integer, intent(in) :: ce
2021-01-25 19:43:17 +05:30
real(pReal), intent(in) :: &
phi
2021-04-08 00:36:29 +05:30
end subroutine homogenization_set_phi
2021-01-25 19:43:17 +05:30
2019-06-15 20:12:16 +05:30
end interface
2019-06-15 20:12:16 +05:30
public :: &
homogenization_init, &
materialpoint_stressAndItsTangent, &
2021-04-09 03:10:20 +05:30
homogenization_mu_T, &
homogenization_K_T, &
2021-04-08 02:11:49 +05:30
homogenization_f_T, &
2021-04-11 12:02:13 +05:30
homogenization_thermal_setfield, &
2021-04-08 17:01:21 +05:30
homogenization_mu_phi, &
2021-04-11 12:02:13 +05:30
homogenization_K_phi, &
2021-04-08 17:01:21 +05:30
homogenization_f_phi, &
2021-04-08 00:36:29 +05:30
homogenization_set_phi, &
2020-12-23 11:28:54 +05:30
homogenization_forward, &
2020-12-30 04:44:48 +05:30
homogenization_results, &
homogenization_restartRead, &
2021-02-12 20:01:43 +05:30
homogenization_restartWrite, &
THERMAL_CONDUCTION_ID, &
DAMAGE_NONLOCAL_ID
2012-08-10 21:28:17 +05:30
contains
!--------------------------------------------------------------------------------------------------
!> @brief module initialization
!--------------------------------------------------------------------------------------------------
2021-01-21 01:24:31 +05:30
subroutine homogenization_init()
2020-06-16 22:17:19 +05:30
class (tNode) , pointer :: &
num_homog, &
num_homogGeneric
2020-07-02 01:50:22 +05:30
2020-12-16 15:51:24 +05:30
print'(/,a)', ' <<<+- homogenization init -+>>>'; flush(IO_STDOUT)
2021-02-12 20:01:43 +05:30
allocate(homogState (size(material_name_homogenization)))
allocate(damageState_h (size(material_name_homogenization)))
call material_parseHomogenization()
num_homog => config_numerics%get('homogenization',defaultVal=emptyDict)
2020-06-16 22:17:19 +05:30
num_homogGeneric => num_homog%get('generic',defaultVal=emptyDict)
num%nMPstate = num_homogGeneric%get_asInt('nMPstate',defaultVal=10)
if (num%nMPstate < 1) call IO_error(301,ext_msg='nMPstate')
2020-12-16 15:51:24 +05:30
2021-02-09 03:51:53 +05:30
call mechanical_init(num_homog)
2021-01-08 02:45:18 +05:30
call thermal_init()
2021-01-21 01:24:31 +05:30
call damage_init()
2020-12-16 15:51:24 +05:30
2012-08-10 21:28:17 +05:30
end subroutine homogenization_init
2012-08-10 21:28:17 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief parallelized calculation of stress and corresponding tangent at material points
!--------------------------------------------------------------------------------------------------
2020-12-27 16:18:02 +05:30
subroutine materialpoint_stressAndItsTangent(dt,FEsolving_execIP,FEsolving_execElem)
2019-06-15 20:12:16 +05:30
real(pReal), intent(in) :: dt !< time increment
2020-12-27 16:18:02 +05:30
integer, dimension(2), intent(in) :: FEsolving_execElem, FEsolving_execIP
2019-06-15 20:12:16 +05:30
integer :: &
NiterationMPstate, &
2020-12-23 18:54:44 +05:30
ip, & !< integration point number
el, & !< element number
2021-04-06 15:08:44 +05:30
myNgrains, co, ce, ho, en, ph
2020-12-23 21:44:26 +05:30
logical :: &
2020-04-15 16:39:05 +05:30
converged
2020-12-23 21:44:26 +05:30
logical, dimension(2) :: &
2020-04-15 16:39:05 +05:30
doneAndHappy
2020-07-03 20:15:11 +05:30
2021-01-17 14:46:56 +05:30
!$OMP PARALLEL
2021-04-06 15:08:44 +05:30
!$OMP DO PRIVATE(ce,en,ho,myNgrains,NiterationMPstate,converged,doneAndHappy)
2020-12-23 18:54:44 +05:30
do el = FEsolving_execElem(1),FEsolving_execElem(2)
2020-12-27 16:18:02 +05:30
ho = material_homogenizationAt(el)
myNgrains = homogenization_Nconstituents(ho)
2020-12-23 21:44:26 +05:30
do ip = FEsolving_execIP(1),FEsolving_execIP(2)
2021-01-24 14:09:40 +05:30
ce = (el-1)*discretization_nIPs + ip
2021-04-06 15:08:44 +05:30
en = material_homogenizationEntry(ce)
2021-01-24 14:09:40 +05:30
2021-02-09 03:51:53 +05:30
call phase_restore(ce,.false.) ! wrong name (is more a forward function)
2020-03-17 05:09:32 +05:30
2021-04-06 15:08:44 +05:30
if(homogState(ho)%sizeState > 0) homogState(ho)%state(:,en) = homogState(ho)%state0(:,en)
if(damageState_h(ho)%sizeState > 0) damageState_h(ho)%state(:,en) = damageState_h(ho)%state0(:,en)
2021-02-12 14:01:02 +05:30
call damage_partition(ce)
2019-06-15 20:12:16 +05:30
2021-01-25 18:19:29 +05:30
doneAndHappy = [.false.,.true.]
2021-01-25 18:19:29 +05:30
NiterationMPstate = 0
convergenceLooping: do while (.not. (terminallyIll .or. doneAndHappy(1)) &
.and. NiterationMPstate < num%nMPstate)
NiterationMPstate = NiterationMPstate + 1
2019-06-15 20:12:16 +05:30
2021-03-18 12:33:56 +05:30
call mechanical_partition(homogenization_F(1:3,1:3,ce),ce)
converged = .true.
do co = 1, myNgrains
converged = converged .and. crystallite_stress(dt,co,ip,el)
enddo
2021-03-18 12:37:10 +05:30
if (converged) then
2021-03-18 12:33:56 +05:30
doneAndHappy = mechanical_updateState(dt,homogenization_F(1:3,1:3,ce),ce)
converged = all(doneAndHappy)
2021-03-18 12:37:10 +05:30
else
doneAndHappy = [.true.,.false.]
2021-01-25 18:19:29 +05:30
endif
2020-03-17 05:09:32 +05:30
2021-01-25 18:19:29 +05:30
enddo convergenceLooping
if (.not. converged) then
if (.not. terminallyIll) print*, ' Integration point ', ip,' at element ', el, ' terminally ill'
terminallyIll = .true.
endif
2020-12-23 19:07:12 +05:30
enddo
2020-12-24 17:14:26 +05:30
enddo
2021-01-17 14:46:56 +05:30
!$OMP END DO
2020-12-23 19:07:12 +05:30
if (.not. terminallyIll) then
2021-01-17 14:46:56 +05:30
!$OMP DO PRIVATE(ho,ph,ce)
do el = FEsolving_execElem(1),FEsolving_execElem(2)
if (terminallyIll) continue
ho = material_homogenizationAt(el)
do ip = FEsolving_execIP(1),FEsolving_execIP(2)
ce = (el-1)*discretization_nIPs + ip
2021-01-24 17:56:01 +05:30
call thermal_partition(ce)
do co = 1, homogenization_Nconstituents(ho)
ph = material_phaseAt(co,el)
if (.not. thermal_stress(dt,ph,material_phaseMemberAt(co,ip,el))) then
2021-04-13 01:08:42 +05:30
if (.not. terminallyIll) & ! so first signals terminally ill...
print*, ' Integration point ', ip,' at element ', el, ' terminally ill'
terminallyIll = .true. ! ...and kills all others
endif
enddo
enddo
enddo
2021-01-17 14:46:56 +05:30
!$OMP END DO
2021-02-26 02:12:40 +05:30
!$OMP DO PRIVATE(ho,ce)
2020-12-23 18:54:44 +05:30
elementLooping3: do el = FEsolving_execElem(1),FEsolving_execElem(2)
2020-12-27 16:18:02 +05:30
ho = material_homogenizationAt(el)
2020-12-23 18:54:44 +05:30
IpLooping3: do ip = FEsolving_execIP(1),FEsolving_execIP(2)
2021-02-26 02:12:40 +05:30
ce = (el-1)*discretization_nIPs + ip
do co = 1, homogenization_Nconstituents(ho)
2020-12-27 16:18:02 +05:30
call crystallite_orientations(co,ip,el)
enddo
2021-02-26 02:12:40 +05:30
call mechanical_homogenize(dt,ce)
2020-12-23 14:35:02 +05:30
enddo IpLooping3
enddo elementLooping3
2021-01-17 14:46:56 +05:30
!$OMP END DO
2019-06-15 20:12:16 +05:30
else
2020-09-18 02:27:56 +05:30
print'(/,a,/)', ' << HOMOG >> Material Point terminally ill'
2019-06-15 20:12:16 +05:30
endif
2021-01-17 14:46:56 +05:30
!$OMP END PARALLEL
2012-08-10 21:28:17 +05:30
end subroutine materialpoint_stressAndItsTangent
2019-04-30 22:15:16 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief writes homogenization results to HDF5 output file
!--------------------------------------------------------------------------------------------------
subroutine homogenization_results
2020-03-17 05:09:32 +05:30
2020-12-28 13:57:48 +05:30
integer :: ho
2020-11-06 05:29:12 +05:30
character(len=:), allocatable :: group_base,group
2020-03-17 05:09:32 +05:30
2020-12-28 13:48:20 +05:30
2020-12-22 13:15:01 +05:30
call results_closeGroup(results_addGroup('current/homogenization/'))
2020-03-17 05:09:32 +05:30
2020-12-28 13:57:48 +05:30
do ho=1,size(material_name_homogenization)
group_base = 'current/homogenization/'//trim(material_name_homogenization(ho))
2019-12-10 11:44:39 +05:30
call results_closeGroup(results_addGroup(group_base))
2020-03-17 05:09:32 +05:30
2021-02-09 03:51:53 +05:30
call mechanical_results(group_base,ho)
2020-12-28 13:57:48 +05:30
select case(damage_type(ho))
2019-12-10 11:44:39 +05:30
case(DAMAGE_NONLOCAL_ID)
2021-03-25 23:52:59 +05:30
group = trim(group_base)//'/damage'
call results_closeGroup(results_addGroup(group))
2021-04-07 10:56:54 +05:30
call damage_results(ho,group)
2019-12-10 11:44:39 +05:30
end select
2020-03-17 05:09:32 +05:30
2020-12-28 13:57:48 +05:30
select case(thermal_type(ho))
2019-12-11 00:55:19 +05:30
case(THERMAL_CONDUCTION_ID)
2021-03-25 23:52:59 +05:30
group = trim(group_base)//'/thermal'
call results_closeGroup(results_addGroup(group))
2021-04-07 10:56:54 +05:30
call thermal_results(ho,group)
2019-12-11 00:55:19 +05:30
end select
2020-03-17 05:09:32 +05:30
2019-12-10 11:44:39 +05:30
enddo
2019-12-19 00:35:51 +05:30
2019-04-30 22:15:16 +05:30
end subroutine homogenization_results
2020-12-23 11:28:54 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Forward data after successful increment.
! ToDo: Any guessing for the current states possible?
!--------------------------------------------------------------------------------------------------
subroutine homogenization_forward
2020-12-23 14:35:02 +05:30
integer :: ho
2020-12-23 11:28:54 +05:30
2020-12-28 13:48:20 +05:30
2020-12-23 11:28:54 +05:30
do ho = 1, size(material_name_homogenization)
homogState (ho)%state0 = homogState (ho)%state
2021-02-12 14:01:02 +05:30
if(damageState_h(ho)%sizeState > 0) &
damageState_h(ho)%state0 = damageState_h(ho)%state
2020-12-23 11:28:54 +05:30
enddo
end subroutine homogenization_forward
2020-12-30 04:44:48 +05:30
!--------------------------------------------------------------------------------------------------
!--------------------------------------------------------------------------------------------------
subroutine homogenization_restartWrite(fileHandle)
integer(HID_T), intent(in) :: fileHandle
integer(HID_T), dimension(2) :: groupHandle
integer :: ho
groupHandle(1) = HDF5_addGroup(fileHandle,'homogenization')
do ho = 1, size(material_name_homogenization)
groupHandle(2) = HDF5_addGroup(groupHandle(1),material_name_homogenization(ho))
2021-02-03 23:10:52 +05:30
call HDF5_write(groupHandle(2),homogState(ho)%state,'omega') ! ToDo: should be done by mech
2020-12-30 04:44:48 +05:30
call HDF5_closeGroup(groupHandle(2))
enddo
call HDF5_closeGroup(groupHandle(1))
end subroutine homogenization_restartWrite
!--------------------------------------------------------------------------------------------------
!--------------------------------------------------------------------------------------------------
subroutine homogenization_restartRead(fileHandle)
integer(HID_T), intent(in) :: fileHandle
integer(HID_T), dimension(2) :: groupHandle
integer :: ho
groupHandle(1) = HDF5_openGroup(fileHandle,'homogenization')
do ho = 1, size(material_name_homogenization)
groupHandle(2) = HDF5_openGroup(groupHandle(1),material_name_homogenization(ho))
2021-02-24 22:43:48 +05:30
call HDF5_read(groupHandle(2),homogState(ho)%state0,'omega') ! ToDo: should be done by mech
2020-12-30 04:44:48 +05:30
call HDF5_closeGroup(groupHandle(2))
enddo
call HDF5_closeGroup(groupHandle(1))
end subroutine homogenization_restartRead
2021-02-12 20:01:43 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief parses the homogenization part from the material configuration
!--------------------------------------------------------------------------------------------------
subroutine material_parseHomogenization
class(tNode), pointer :: &
material_homogenization, &
homog, &
homogThermal, &
homogDamage
integer :: h
material_homogenization => config_material%get('homogenization')
allocate(thermal_type(size(material_name_homogenization)), source=THERMAL_isothermal_ID)
allocate(damage_type (size(material_name_homogenization)), source=DAMAGE_none_ID)
2021-02-12 20:01:43 +05:30
do h=1, size(material_name_homogenization)
homog => material_homogenization%get(h)
if (homog%contains('thermal')) then
2021-02-12 20:01:43 +05:30
homogThermal => homog%get('thermal')
select case (homogThermal%get_asString('type'))
case('pass')
2021-02-12 20:01:43 +05:30
thermal_type(h) = THERMAL_conduction_ID
case default
call IO_error(500,ext_msg=homogThermal%get_asString('type'))
end select
endif
if (homog%contains('damage')) then
2021-02-12 20:01:43 +05:30
homogDamage => homog%get('damage')
select case (homogDamage%get_asString('type'))
2021-03-01 10:46:16 +05:30
case('pass')
2021-02-12 20:01:43 +05:30
damage_type(h) = DAMAGE_nonlocal_ID
case default
call IO_error(500,ext_msg=homogDamage%get_asString('type'))
end select
endif
enddo
end subroutine material_parseHomogenization
2012-08-10 21:28:17 +05:30
end module homogenization