DAMASK_EICMD/src/homogenization.f90

446 lines
21 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 constitutive
use FEsolving
use discretization
use thermal_isothermal
use thermal_conduction
use damage_none
use damage_nonlocal
use results
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
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
real(pReal) :: &
subStepMinHomog, & !< minimum (relative) size of sub-step allowed during cutback in homogenization
subStepSizeHomog, & !< size of first substep when cutback in homogenization
stepIncreaseHomog !< increase of next substep size when previous substep converged in homogenization
end type tNumerics
type(tNumerics) :: num
2020-07-02 01:50:22 +05:30
type :: tDebugOptions
logical :: &
basic, &
extensive, &
selective
integer :: &
element, &
ip, &
grain
end type tDebugOptions
type(tDebugOptions) :: debugHomog
2020-07-02 01:50:22 +05:30
2020-12-16 15:51:24 +05:30
!--------------------------------------------------------------------------------------------------
interface
2020-03-17 05:09:32 +05:30
2020-12-16 15:51:24 +05:30
module subroutine mech_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
end subroutine mech_init
module subroutine mech_partition(subF,ip,el)
real(pReal), intent(in), dimension(3,3) :: &
subF
integer, intent(in) :: &
ip, & !< integration point
el !< element number
end subroutine mech_partition
2020-03-17 05:09:32 +05:30
2020-12-16 15:51:24 +05:30
module subroutine mech_homogenize(ip,el)
integer, intent(in) :: &
ip, & !< integration point
el !< element number
end subroutine mech_homogenize
2020-03-17 05:09:32 +05:30
2020-12-16 15:51:24 +05:30
module subroutine mech_results(group_base,h)
2020-03-17 05:09:32 +05:30
2020-12-16 15:51:24 +05:30
character(len=*), intent(in) :: group_base
integer, intent(in) :: h
2020-03-17 05:09:32 +05:30
2020-12-16 15:51:24 +05:30
end subroutine mech_results
2020-03-17 05:09:32 +05:30
2020-12-16 15:51:24 +05:30
! -------- ToDo ---------------------------------------------------------
2020-07-02 01:50:22 +05:30
module function mech_RGC_updateState(P,F,F0,avgF,dt,dPdF,ip,el)
2019-06-15 20:12:16 +05:30
logical, dimension(2) :: mech_RGC_updateState
2020-03-17 05:09:32 +05:30
real(pReal), dimension(:,:,:), intent(in) :: &
2019-06-15 20:12:16 +05:30
P,& !< partitioned stresses
F,& !< partitioned deformation gradients
F0 !< partitioned initial deformation gradients
2020-06-26 23:42:05 +05:30
real(pReal), dimension(:,:,:,:,:), intent(in) :: dPdF !< partitioned stiffnesses
real(pReal), dimension(3,3), intent(in) :: avgF !< average F
real(pReal), intent(in) :: dt !< time increment
integer, intent(in) :: &
2020-06-26 23:42:05 +05:30
ip, & !< integration point number
el !< element number
2019-06-15 20:12:16 +05:30
end function mech_RGC_updateState
end interface
2020-12-16 15:51:24 +05:30
! -----------------------------------------------------------------------
2019-06-15 20:12:16 +05:30
public :: &
homogenization_init, &
materialpoint_stressAndItsTangent, &
2020-12-23 11:28:54 +05:30
homogenization_forward, &
2019-06-15 20:12:16 +05:30
homogenization_results
2012-08-10 21:28:17 +05:30
contains
!--------------------------------------------------------------------------------------------------
!> @brief module initialization
!--------------------------------------------------------------------------------------------------
subroutine homogenization_init
2020-06-16 22:17:19 +05:30
class (tNode) , pointer :: &
num_homog, &
num_homogGeneric, &
debug_homogenization
2020-07-02 01:50:22 +05:30
2020-12-16 15:51:24 +05:30
print'(/,a)', ' <<<+- homogenization init -+>>>'; flush(IO_STDOUT)
debug_homogenization => config_debug%get('homogenization', defaultVal=emptyList)
2020-08-13 00:44:00 +05:30
debugHomog%basic = debug_homogenization%contains('basic')
debugHomog%extensive = debug_homogenization%contains('extensive')
debugHomog%selective = debug_homogenization%contains('selective')
debugHomog%element = config_debug%get_asInt('element',defaultVal = 1)
debugHomog%ip = config_debug%get_asInt('integrationpoint',defaultVal = 1)
debugHomog%grain = config_debug%get_asInt('grain',defaultVal = 1)
2020-07-02 01:50:22 +05:30
if (debugHomog%grain < 1 &
.or. debugHomog%grain > homogenization_Nconstituents(material_homogenizationAt(debugHomog%element))) &
call IO_error(602,ext_msg='constituent', el=debugHomog%element, g=debugHomog%grain)
2020-07-02 01:50:22 +05:30
num_homog => config_numerics%get('homogenization',defaultVal=emptyDict)
2020-06-16 22:17:19 +05:30
num_homogGeneric => num_homog%get('generic',defaultVal=emptyDict)
2020-06-26 15:52:33 +05:30
num%nMPstate = num_homogGeneric%get_asInt ('nMPstate', defaultVal=10)
2020-06-16 22:17:19 +05:30
num%subStepMinHomog = num_homogGeneric%get_asFloat('subStepMin', defaultVal=1.0e-3_pReal)
num%subStepSizeHomog = num_homogGeneric%get_asFloat('subStepSize', defaultVal=0.25_pReal)
num%stepIncreaseHomog = num_homogGeneric%get_asFloat('stepIncrease', defaultVal=1.5_pReal)
2020-03-17 05:09:32 +05:30
if (num%nMPstate < 1) call IO_error(301,ext_msg='nMPstate')
if (num%subStepMinHomog <= 0.0_pReal) call IO_error(301,ext_msg='subStepMinHomog')
if (num%subStepSizeHomog <= 0.0_pReal) call IO_error(301,ext_msg='subStepSizeHomog')
if (num%stepIncreaseHomog <= 0.0_pReal) call IO_error(301,ext_msg='stepIncreaseHomog')
2020-12-16 15:51:24 +05:30
call mech_init(num_homog)
if (any(thermal_type == THERMAL_isothermal_ID)) call thermal_isothermal_init
if (any(thermal_type == THERMAL_conduction_ID)) call thermal_conduction_init
if (any(damage_type == DAMAGE_none_ID)) call damage_none_init
if (any(damage_type == DAMAGE_nonlocal_ID)) call damage_nonlocal_init
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
!--------------------------------------------------------------------------------------------------
subroutine materialpoint_stressAndItsTangent(dt)
2019-06-15 20:12:16 +05:30
real(pReal), intent(in) :: dt !< time increment
integer :: &
NiterationHomog, &
NiterationMPstate, &
i, & !< integration point number
e, & !< element number
2020-07-02 01:50:22 +05:30
myNgrains
real(pReal), dimension(discretization_nIPs,discretization_Nelems) :: &
2020-04-15 16:39:05 +05:30
subFrac, &
subStep
logical, dimension(discretization_nIPs,discretization_Nelems) :: &
2020-04-15 16:39:05 +05:30
requested, &
converged
logical, dimension(2,discretization_nIPs,discretization_Nelems) :: &
2020-04-15 16:39:05 +05:30
doneAndHappy
integer :: m
2020-07-03 20:15:11 +05:30
!--------------------------------------------------------------------------------------------------
2020-04-15 12:23:25 +05:30
! initialize restoration points
2019-06-15 20:12:16 +05:30
do e = FEsolving_execElem(1),FEsolving_execElem(2)
do i = FEsolving_execIP(1),FEsolving_execIP(2);
2019-06-15 20:12:16 +05:30
2020-12-22 04:03:32 +05:30
call constitutive_initializeRestorationPoints(i,e)
2019-06-15 20:12:16 +05:30
2020-04-15 16:39:05 +05:30
subFrac(i,e) = 0.0_pReal
converged(i,e) = .false. ! pretend failed step ...
subStep(i,e) = 1.0_pReal/num%subStepSizeHomog ! ... larger then the requested calculation
requested(i,e) = .true. ! everybody requires calculation
2020-03-17 05:09:32 +05:30
2019-06-15 20:12:16 +05:30
if (homogState(material_homogenizationAt(e))%sizeState > 0) &
homogState(material_homogenizationAt(e))%subState0(:,material_homogenizationMemberAt(i,e)) = &
2020-04-15 12:23:25 +05:30
homogState(material_homogenizationAt(e))%State0( :,material_homogenizationMemberAt(i,e))
2019-06-15 20:12:16 +05:30
if (damageState(material_homogenizationAt(e))%sizeState > 0) &
damageState(material_homogenizationAt(e))%subState0(:,material_homogenizationMemberAt(i,e)) = &
2020-04-15 12:23:25 +05:30
damageState(material_homogenizationAt(e))%State0( :,material_homogenizationMemberAt(i,e))
2019-06-15 20:12:16 +05:30
enddo
enddo
2020-03-17 05:09:32 +05:30
2019-06-15 20:12:16 +05:30
NiterationHomog = 0
2019-06-15 20:12:16 +05:30
cutBackLooping: do while (.not. terminallyIll .and. &
any(subStep(FEsolving_execIP(1):FEsolving_execIP(2),&
FEsolving_execElem(1):FEsolving_execElem(2)) > num%subStepMinHomog))
!$OMP PARALLEL DO PRIVATE(m)
2019-06-15 20:12:16 +05:30
elementLooping1: do e = FEsolving_execElem(1),FEsolving_execElem(2)
myNgrains = homogenization_Nconstituents(material_homogenizationAt(e))
IpLooping1: do i = FEsolving_execIP(1),FEsolving_execIP(2)
2020-04-15 16:39:05 +05:30
if (converged(i,e)) then
subFrac(i,e) = subFrac(i,e) + subStep(i,e)
subStep(i,e) = min(1.0_pReal-subFrac(i,e),num%stepIncreaseHomog*subStep(i,e)) ! introduce flexibility for step increase/acceleration
2019-06-15 20:12:16 +05:30
2020-04-15 16:39:05 +05:30
steppingNeeded: if (subStep(i,e) > num%subStepMinHomog) then
2019-06-15 20:12:16 +05:30
2020-04-15 12:23:25 +05:30
! wind forward grain starting point
2020-12-22 14:53:46 +05:30
call constitutive_windForward(i,e)
2019-06-15 20:12:16 +05:30
if(homogState(material_homogenizationAt(e))%sizeState > 0) &
homogState(material_homogenizationAt(e))%subState0(:,material_homogenizationMemberAt(i,e)) = &
homogState(material_homogenizationAt(e))%State (:,material_homogenizationMemberAt(i,e))
2019-06-15 20:12:16 +05:30
if(damageState(material_homogenizationAt(e))%sizeState > 0) &
damageState(material_homogenizationAt(e))%subState0(:,material_homogenizationMemberAt(i,e)) = &
damageState(material_homogenizationAt(e))%State (:,material_homogenizationMemberAt(i,e))
2020-03-17 05:09:32 +05:30
2019-06-15 20:12:16 +05:30
endif steppingNeeded
2020-04-15 16:39:05 +05:30
else
if ( (myNgrains == 1 .and. subStep(i,e) <= 1.0 ) .or. & ! single grain already tried internal subStepping in crystallite
num%subStepSizeHomog * subStep(i,e) <= num%subStepMinHomog ) then ! would require too small subStep
! cutback makes no sense
2019-06-15 20:12:16 +05:30
if (.not. terminallyIll) then ! so first signals terminally ill...
2020-09-18 02:27:56 +05:30
print*, ' Integration point ', i,' at element ', e, ' terminally ill'
2019-06-15 20:12:16 +05:30
endif
terminallyIll = .true. ! ...and kills all others
2019-06-15 20:12:16 +05:30
else ! cutback makes sense
2020-04-15 16:39:05 +05:30
subStep(i,e) = num%subStepSizeHomog * subStep(i,e) ! crystallite had severe trouble, so do a significant cutback
2020-09-28 21:26:48 +05:30
call crystallite_restore(i,e,subStep(i,e) < 1.0_pReal)
2020-12-20 13:57:37 +05:30
call constitutive_restore(i,e)
2019-06-15 20:12:16 +05:30
if(homogState(material_homogenizationAt(e))%sizeState > 0) &
homogState(material_homogenizationAt(e))%State( :,material_homogenizationMemberAt(i,e)) = &
homogState(material_homogenizationAt(e))%subState0(:,material_homogenizationMemberAt(i,e))
2019-06-15 20:12:16 +05:30
if(damageState(material_homogenizationAt(e))%sizeState > 0) &
damageState(material_homogenizationAt(e))%State( :,material_homogenizationMemberAt(i,e)) = &
damageState(material_homogenizationAt(e))%subState0(:,material_homogenizationMemberAt(i,e))
2019-06-15 20:12:16 +05:30
endif
2020-04-15 16:39:05 +05:30
endif
2019-06-15 20:12:16 +05:30
2020-04-15 16:39:05 +05:30
if (subStep(i,e) > num%subStepMinHomog) then
requested(i,e) = .true.
doneAndHappy(1:2,i,e) = [.false.,.true.]
2019-06-15 20:12:16 +05:30
endif
enddo IpLooping1
enddo elementLooping1
!$OMP END PARALLEL DO
NiterationMPstate = 0
convergenceLooping: do while (.not. terminallyIll .and. &
2020-04-15 16:39:05 +05:30
any( requested(:,FEsolving_execELem(1):FEsolving_execElem(2)) &
.and. .not. doneAndHappy(1,:,FEsolving_execELem(1):FEsolving_execElem(2)) &
2019-06-15 20:12:16 +05:30
) .and. &
2020-03-17 05:09:32 +05:30
NiterationMPstate < num%nMPstate)
2019-06-15 20:12:16 +05:30
NiterationMPstate = NiterationMPstate + 1
!--------------------------------------------------------------------------------------------------
! deformation partitioning
!$OMP PARALLEL DO PRIVATE(myNgrains,m)
2019-06-15 20:12:16 +05:30
elementLooping2: do e = FEsolving_execElem(1),FEsolving_execElem(2)
myNgrains = homogenization_Nconstituents(material_homogenizationAt(e))
IpLooping2: do i = FEsolving_execIP(1),FEsolving_execIP(2)
2020-04-15 16:39:05 +05:30
if(requested(i,e) .and. .not. doneAndHappy(1,i,e)) then ! requested but not yet done
m = (e-1)*discretization_nIPs + i
call mech_partition(homogenization_F0(1:3,1:3,m) &
+ (homogenization_F(1:3,1:3,m)-homogenization_F0(1:3,1:3,m))&
2020-04-16 01:26:20 +05:30
*(subStep(i,e)+subFrac(i,e)), &
i,e)
2020-04-15 16:39:05 +05:30
crystallite_dt(1:myNgrains,i,e) = dt*subStep(i,e) ! propagate materialpoint dt to grains
2019-06-15 20:12:16 +05:30
crystallite_requested(1:myNgrains,i,e) = .true. ! request calculation for constituents
else
crystallite_requested(1:myNgrains,i,e) = .false. ! calculation for constituents not required anymore
endif
enddo IpLooping2
enddo elementLooping2
!$OMP END PARALLEL DO
!--------------------------------------------------------------------------------------------------
! crystallite integration
2020-04-15 16:39:05 +05:30
converged = crystallite_stress() !ToDo: MD not sure if that is the best logic
!--------------------------------------------------------------------------------------------------
! state update
!$OMP PARALLEL DO PRIVATE(m)
2019-06-15 20:12:16 +05:30
elementLooping3: do e = FEsolving_execElem(1),FEsolving_execElem(2)
IpLooping3: do i = FEsolving_execIP(1),FEsolving_execIP(2)
2020-04-15 16:39:05 +05:30
if (requested(i,e) .and. .not. doneAndHappy(1,i,e)) then
if (.not. converged(i,e)) then
doneAndHappy(1:2,i,e) = [.true.,.false.]
2019-06-15 20:12:16 +05:30
else
m = (e-1)*discretization_nIPs + i
2020-04-16 01:26:20 +05:30
doneAndHappy(1:2,i,e) = updateState(dt*subStep(i,e), &
homogenization_F0(1:3,1:3,m) &
+ (homogenization_F(1:3,1:3,m)-homogenization_F0(1:3,1:3,m)) &
2020-04-16 01:26:20 +05:30
*(subStep(i,e)+subFrac(i,e)), &
i,e)
2020-04-15 16:39:05 +05:30
converged(i,e) = all(doneAndHappy(1:2,i,e)) ! converged if done and happy
2019-06-15 20:12:16 +05:30
endif
endif
enddo IpLooping3
enddo elementLooping3
!$OMP END PARALLEL DO
2020-03-17 05:09:32 +05:30
2019-06-15 20:12:16 +05:30
enddo convergenceLooping
2020-03-17 05:09:32 +05:30
2019-06-15 20:12:16 +05:30
NiterationHomog = NiterationHomog + 1
2020-03-17 05:09:32 +05:30
2019-06-15 20:12:16 +05:30
enddo cutBackLooping
2020-03-17 05:09:32 +05:30
2019-06-15 20:12:16 +05:30
if (.not. terminallyIll ) then
call crystallite_orientations() ! calculate crystal orientations
!$OMP PARALLEL DO
elementLooping4: do e = FEsolving_execElem(1),FEsolving_execElem(2)
IpLooping4: do i = FEsolving_execIP(1),FEsolving_execIP(2)
2020-12-16 15:51:24 +05:30
call mech_homogenize(i,e)
2019-06-15 20:12:16 +05:30
enddo IpLooping4
enddo elementLooping4
!$OMP END PARALLEL DO
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
2012-08-10 21:28:17 +05:30
end subroutine materialpoint_stressAndItsTangent
2012-08-10 21:28:17 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief update the internal state of the homogenization scheme and tell whether "done" and
2012-08-10 21:28:17 +05:30
!> "happy" with result
!--------------------------------------------------------------------------------------------------
2020-04-14 11:15:39 +05:30
function updateState(subdt,subF,ip,el)
2020-04-14 13:19:03 +05:30
real(pReal), intent(in) :: &
subdt !< current time step
real(pReal), intent(in), dimension(3,3) :: &
subF
integer, intent(in) :: &
ip, & !< integration point
el !< element number
integer :: c
2020-04-14 13:19:03 +05:30
logical, dimension(2) :: updateState
real(pReal) :: dPdFs(3,3,3,3,homogenization_Nconstituents(material_homogenizationAt(el)))
2020-07-03 20:15:11 +05:30
2020-04-14 13:19:03 +05:30
updateState = .true.
chosenHomogenization: select case(homogenization_type(material_homogenizationAt(el)))
case (HOMOGENIZATION_RGC_ID) chosenHomogenization
do c=1,homogenization_Nconstituents(material_homogenizationAt(el))
dPdFs(:,:,:,:,c) = crystallite_stressTangent(c,ip,el)
enddo
2020-04-14 13:19:03 +05:30
updateState = &
updateState .and. &
mech_RGC_updateState(crystallite_P(1:3,1:3,1:homogenization_Nconstituents(material_homogenizationAt(el)),ip,el), &
crystallite_partitionedF(1:3,1:3,1:homogenization_Nconstituents(material_homogenizationAt(el)),ip,el), &
crystallite_partitionedF0(1:3,1:3,1:homogenization_Nconstituents(material_homogenizationAt(el)),ip,el),&
2020-04-14 13:19:03 +05:30
subF,&
subdt, &
dPdFs, &
2020-04-14 13:19:03 +05:30
ip, &
2020-07-02 01:50:22 +05:30
el)
2020-04-14 13:19:03 +05:30
end select chosenHomogenization
2019-01-13 14:18:47 +05:30
end function updateState
2019-04-30 22:15:16 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief writes homogenization results to HDF5 output file
!--------------------------------------------------------------------------------------------------
subroutine homogenization_results
use material, only: &
material_homogenization_type => homogenization_type
2020-03-17 05:09:32 +05:30
2019-04-30 22:15:16 +05:30
integer :: p
2020-11-06 05:29:12 +05:30
character(len=:), allocatable :: group_base,group
2020-03-17 05:09:32 +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-08-15 19:32:10 +05:30
do p=1,size(material_name_homogenization)
group_base = 'current/homogenization/'//trim(material_name_homogenization(p))
2019-12-10 11:44:39 +05:30
call results_closeGroup(results_addGroup(group_base))
2020-03-17 05:09:32 +05:30
2020-12-16 15:51:24 +05:30
call mech_results(group_base,p)
2019-12-10 11:44:39 +05:30
group = trim(group_base)//'/damage'
call results_closeGroup(results_addGroup(group))
select case(damage_type(p))
case(DAMAGE_NONLOCAL_ID)
call damage_nonlocal_results(p,group)
end select
2020-03-17 05:09:32 +05:30
2019-12-10 11:44:39 +05:30
group = trim(group_base)//'/thermal'
call results_closeGroup(results_addGroup(group))
2019-12-11 00:55:19 +05:30
select case(thermal_type(p))
case(THERMAL_CONDUCTION_ID)
call thermal_conduction_results(p,group)
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
integer :: ho
do ho = 1, size(material_name_homogenization)
homogState (ho)%state0 = homogState (ho)%state
damageState(ho)%state0 = damageState(ho)%state
enddo
end subroutine homogenization_forward
2012-08-10 21:28:17 +05:30
end module homogenization