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
2015-10-15 00:06:19 +05:30
!> @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 crystallite
use FEsolving
use discretization
use thermal_isothermal
use thermal_adiabatic
use thermal_conduction
use damage_none
use damage_local
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 :: &
2019-06-15 20:12:16 +05:30
materialpoint_F0 , & !< def grad of IP at start of FE increment
2020-04-14 13:08:32 +05:30
materialpoint_F !< def grad of IP to be reached at end of FE increment
real ( pReal ) , dimension ( : , : , : , : ) , allocatable , public , protected :: &
2019-06-15 20:12:16 +05:30
materialpoint_P !< first P--K stress of IP
2020-04-14 13:08:32 +05:30
real ( pReal ) , dimension ( : , : , : , : , : , : ) , allocatable , public , protected :: &
2019-06-15 20:12:16 +05:30
materialpoint_dPdF !< tangent of first P--K stress at IP
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
2020-07-02 04:55:24 +05:30
type ( tDebugOptions ) :: debugHomog
2020-07-02 01:50:22 +05:30
2019-06-15 20:12:16 +05:30
interface
2019-04-06 01:12:07 +05:30
2020-07-02 01:50:22 +05:30
module subroutine mech_none_init
2019-06-15 20:12:16 +05:30
end subroutine mech_none_init
2020-03-17 05:09:32 +05:30
2020-07-02 01:50:22 +05:30
module subroutine mech_isostrain_init
2019-06-15 20:12:16 +05:30
end subroutine mech_isostrain_init
2020-03-17 05:09:32 +05:30
2020-07-02 01:50:22 +05:30
module subroutine mech_RGC_init ( num_homogMech )
2020-06-16 22:17:19 +05:30
class ( tNode ) , pointer , intent ( in ) :: &
2020-07-02 01:50:22 +05:30
num_homogMech !< pointer to mechanical homogenization numerics data
2019-06-15 20:12:16 +05:30
end subroutine mech_RGC_init
2020-03-17 05:09:32 +05:30
2019-06-15 20:12:16 +05:30
module subroutine mech_isostrain_partitionDeformation ( F , avgF )
real ( pReal ) , dimension ( : , : , : ) , intent ( out ) :: F !< partitioned deformation gradient
real ( pReal ) , dimension ( 3 , 3 ) , intent ( in ) :: avgF !< average deformation gradient at material point
end subroutine mech_isostrain_partitionDeformation
2020-03-17 05:09:32 +05:30
2020-07-02 01:50:22 +05:30
module subroutine mech_RGC_partitionDeformation ( F , avgF , instance , of )
2019-06-15 20:12:16 +05:30
real ( pReal ) , dimension ( : , : , : ) , intent ( out ) :: F !< partitioned deformation gradient
real ( pReal ) , dimension ( 3 , 3 ) , intent ( in ) :: avgF !< average deformation gradient at material point
integer , intent ( in ) :: &
instance , &
of
end subroutine mech_RGC_partitionDeformation
2020-03-17 05:09:32 +05:30
2019-06-15 20:12:16 +05:30
module subroutine mech_isostrain_averageStressAndItsTangent ( avgP , dAvgPdAvgF , P , dPdF , instance )
real ( pReal ) , dimension ( 3 , 3 ) , intent ( out ) :: avgP !< average stress at material point
real ( pReal ) , dimension ( 3 , 3 , 3 , 3 ) , intent ( out ) :: dAvgPdAvgF !< average stiffness at material point
2020-03-17 05:09:32 +05:30
2019-06-15 20:12:16 +05:30
real ( pReal ) , dimension ( : , : , : ) , intent ( in ) :: P !< partitioned stresses
real ( pReal ) , dimension ( : , : , : , : , : ) , intent ( in ) :: dPdF !< partitioned stiffnesses
2020-03-17 05:09:32 +05:30
integer , intent ( in ) :: instance
2019-06-15 20:12:16 +05:30
end subroutine mech_isostrain_averageStressAndItsTangent
2020-03-17 05:09:32 +05:30
2019-06-15 20:12:16 +05:30
module subroutine mech_RGC_averageStressAndItsTangent ( avgP , dAvgPdAvgF , P , dPdF , instance )
real ( pReal ) , dimension ( 3 , 3 ) , intent ( out ) :: avgP !< average stress at material point
real ( pReal ) , dimension ( 3 , 3 , 3 , 3 ) , intent ( out ) :: dAvgPdAvgF !< average stiffness at material point
2020-03-17 05:09:32 +05:30
2019-06-15 20:12:16 +05:30
real ( pReal ) , dimension ( : , : , : ) , intent ( in ) :: P !< partitioned stresses
real ( pReal ) , dimension ( : , : , : , : , : ) , intent ( in ) :: dPdF !< partitioned stiffnesses
2020-03-17 05:09:32 +05:30
integer , intent ( in ) :: instance
2019-06-15 20:12:16 +05:30
end subroutine mech_RGC_averageStressAndItsTangent
2020-03-17 05:09:32 +05:30
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
2020-06-18 19:36:11 +05:30
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
module subroutine mech_RGC_results ( instance , group )
integer , intent ( in ) :: instance !< homogenization instance
character ( len = * ) , intent ( in ) :: group !< group name in HDF5 file
end subroutine mech_RGC_results
2020-03-17 05:09:32 +05:30
2019-06-15 20:12:16 +05:30
end interface
2013-01-29 15:58:01 +05:30
2019-06-15 20:12:16 +05:30
public :: &
homogenization_init , &
materialpoint_stressAndItsTangent , &
homogenization_results
2012-08-10 21:28:17 +05:30
contains
!--------------------------------------------------------------------------------------------------
!> @brief module initialization
!--------------------------------------------------------------------------------------------------
2015-07-24 20:27:29 +05:30
subroutine homogenization_init
2012-05-08 20:27:06 +05:30
2020-06-16 22:17:19 +05:30
class ( tNode ) , pointer :: &
num_homog , &
num_homogMech , &
2020-06-18 19:36:11 +05:30
num_homogGeneric , &
debug_homogenization
2020-07-02 01:50:22 +05:30
debug_homogenization = > debug_root % get ( 'homogenization' , defaultVal = emptyList )
2020-08-13 00:44:00 +05:30
debugHomog % basic = debug_homogenization % contains ( 'basic' )
debugHomog % extensive = debug_homogenization % contains ( 'extensive' )
2020-07-02 04:55:24 +05:30
debugHomog % selective = debug_homogenization % contains ( 'selective' )
2020-08-13 00:44:00 +05:30
debugHomog % element = debug_root % get_asInt ( 'element' , defaultVal = 1 )
debugHomog % ip = debug_root % get_asInt ( 'integrationpoint' , defaultVal = 1 )
debugHomog % grain = debug_root % get_asInt ( 'grain' , defaultVal = 1 )
2020-07-02 01:50:22 +05:30
2020-07-02 04:55:24 +05:30
if ( debugHomog % grain < 1 &
. or . debugHomog % grain > homogenization_Ngrains ( 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
2020-06-26 23:42:05 +05:30
num_homog = > numerics_root % get ( 'homogenization' , defaultVal = emptyDict )
num_homogMech = > num_homog % get ( 'mech' , defaultVal = emptyDict )
2020-06-16 22:17:19 +05:30
num_homogGeneric = > num_homog % get ( 'generic' , defaultVal = emptyDict )
2020-07-02 01:50:22 +05:30
if ( any ( homogenization_type == HOMOGENIZATION_NONE_ID ) ) call mech_none_init
if ( any ( homogenization_type == HOMOGENIZATION_ISOSTRAIN_ID ) ) call mech_isostrain_init
if ( any ( homogenization_type == HOMOGENIZATION_RGC_ID ) ) call mech_RGC_init ( num_homogMech )
2015-09-24 23:20:11 +05:30
2019-06-15 20:12:16 +05:30
if ( any ( thermal_type == THERMAL_isothermal_ID ) ) call thermal_isothermal_init
if ( any ( thermal_type == THERMAL_adiabatic_ID ) ) call thermal_adiabatic_init
if ( any ( thermal_type == THERMAL_conduction_ID ) ) call thermal_conduction_init
2015-09-24 23:20:11 +05:30
2019-06-15 20:12:16 +05:30
if ( any ( damage_type == DAMAGE_none_ID ) ) call damage_none_init
if ( any ( damage_type == DAMAGE_local_ID ) ) call damage_local_init
if ( any ( damage_type == DAMAGE_nonlocal_ID ) ) call damage_nonlocal_init
2015-09-24 23:20:11 +05:30
2018-06-27 00:20:06 +05:30
2012-08-10 21:28:17 +05:30
!--------------------------------------------------------------------------------------------------
! allocate and initialize global variables
2019-06-15 20:12:16 +05:30
allocate ( materialpoint_dPdF ( 3 , 3 , 3 , 3 , discretization_nIP , discretization_nElem ) , source = 0.0_pReal )
materialpoint_F0 = spread ( spread ( math_I3 , 3 , discretization_nIP ) , 4 , discretization_nElem ) ! initialize to identity
materialpoint_F = materialpoint_F0 ! initialize to identity
allocate ( materialpoint_P ( 3 , 3 , discretization_nIP , discretization_nElem ) , source = 0.0_pReal )
2013-12-16 17:28:03 +05:30
2020-01-14 01:22:58 +05:30
write ( 6 , '(/,a)' ) ' <<<+- homogenization init -+>>>' ; flush ( 6 )
2015-10-15 00:06:19 +05:30
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' )
2012-08-10 21:28:17 +05:30
end subroutine homogenization_init
2009-05-07 21:57:36 +05:30
2012-08-10 21:28:17 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief parallelized calculation of stress and corresponding tangent at material points
!--------------------------------------------------------------------------------------------------
subroutine materialpoint_stressAndItsTangent ( updateJaco , dt )
2015-10-15 00:06:19 +05:30
2019-06-15 20:12:16 +05:30
real ( pReal ) , intent ( in ) :: dt !< time increment
logical , intent ( in ) :: updateJaco !< initiating Jacobian update
integer :: &
NiterationHomog , &
NiterationMPstate , &
g , & !< grain number
i , & !< integration point number
e , & !< element number
mySource , &
2020-07-02 01:50:22 +05:30
myNgrains
2020-04-14 13:13:18 +05:30
real ( pReal ) , dimension ( discretization_nIP , discretization_nElem ) :: &
2020-04-15 16:39:05 +05:30
subFrac , &
subStep
2020-04-14 13:13:18 +05:30
logical , dimension ( discretization_nIP , discretization_nElem ) :: &
2020-04-15 16:39:05 +05:30
requested , &
converged
2020-04-14 13:13:18 +05:30
logical , dimension ( 2 , discretization_nIP , discretization_nElem ) :: &
2020-04-15 16:39:05 +05:30
doneAndHappy
2020-07-03 20:15:11 +05:30
2018-09-20 10:57:12 +05:30
#ifdef DEBUG
2020-06-18 19:36:11 +05:30
2020-07-02 04:55:24 +05:30
if ( debugHomog % basic ) then
write ( 6 , '(/a,i5,1x,i2)' ) '<< HOMOG >> Material Point start at el ip ' , debugHomog % element , debugHomog % ip
2019-06-15 20:12:16 +05:30
write ( 6 , '(a,/,3(12x,3(f14.9,1x)/))' ) '<< HOMOG >> F0' , &
2020-07-02 04:55:24 +05:30
transpose ( materialpoint_F0 ( 1 : 3 , 1 : 3 , debugHomog % ip , debugHomog % element ) )
2019-06-15 20:12:16 +05:30
write ( 6 , '(a,/,3(12x,3(f14.9,1x)/))' ) '<< HOMOG >> F' , &
2020-07-02 04:55:24 +05:30
transpose ( materialpoint_F ( 1 : 3 , 1 : 3 , debugHomog % ip , debugHomog % element ) )
2019-06-15 20:12:16 +05:30
endif
2018-09-20 10:57:12 +05:30
#endif
2009-05-07 21:57:36 +05:30
2013-01-29 15:58:01 +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 )
myNgrains = homogenization_Ngrains ( material_homogenizationAt ( e ) )
2020-01-25 13:54:42 +05:30
do i = FEsolving_execIP ( 1 ) , FEsolving_execIP ( 2 ) ;
2019-06-15 20:12:16 +05:30
do g = 1 , myNgrains
plasticState ( material_phaseAt ( g , e ) ) % partionedState0 ( : , material_phasememberAt ( g , i , e ) ) = &
plasticState ( material_phaseAt ( g , e ) ) % state0 ( : , material_phasememberAt ( g , i , e ) )
do mySource = 1 , phase_Nsources ( material_phaseAt ( g , e ) )
sourceState ( material_phaseAt ( g , e ) ) % p ( mySource ) % partionedState0 ( : , material_phasememberAt ( g , i , e ) ) = &
sourceState ( material_phaseAt ( g , e ) ) % p ( mySource ) % state0 ( : , material_phasememberAt ( g , i , e ) )
enddo
crystallite_partionedFp0 ( 1 : 3 , 1 : 3 , g , i , e ) = crystallite_Fp0 ( 1 : 3 , 1 : 3 , g , i , e )
crystallite_partionedLp0 ( 1 : 3 , 1 : 3 , g , i , e ) = crystallite_Lp0 ( 1 : 3 , 1 : 3 , g , i , e )
crystallite_partionedFi0 ( 1 : 3 , 1 : 3 , g , i , e ) = crystallite_Fi0 ( 1 : 3 , 1 : 3 , g , i , e )
crystallite_partionedLi0 ( 1 : 3 , 1 : 3 , g , i , e ) = crystallite_Li0 ( 1 : 3 , 1 : 3 , g , i , e )
crystallite_partionedF0 ( 1 : 3 , 1 : 3 , g , i , e ) = crystallite_F0 ( 1 : 3 , 1 : 3 , g , i , e )
crystallite_partionedS0 ( 1 : 3 , 1 : 3 , g , i , e ) = crystallite_S0 ( 1 : 3 , 1 : 3 , g , i , e )
enddo
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 ) &
2020-01-23 17:18:18 +05:30
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 ( thermalState ( material_homogenizationAt ( e ) ) % sizeState > 0 ) &
2020-01-23 17:18:18 +05:30
thermalState ( material_homogenizationAt ( e ) ) % subState0 ( : , material_homogenizationMemberAt ( i , e ) ) = &
2020-04-15 12:23:25 +05:30
thermalState ( material_homogenizationAt ( e ) ) % State0 ( : , material_homogenizationMemberAt ( i , e ) )
2020-03-17 05:09:32 +05:30
2019-06-15 20:12:16 +05:30
if ( damageState ( material_homogenizationAt ( e ) ) % sizeState > 0 ) &
2020-01-23 17:18:18 +05:30
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
2015-10-15 00:06:19 +05:30
2019-06-15 20:12:16 +05:30
cutBackLooping : do while ( . not . terminallyIll . and . &
2020-06-16 02:42:49 +05:30
any ( subStep ( FEsolving_execIP ( 1 ) : FEsolving_execIP ( 2 ) , &
FEsolving_execElem ( 1 ) : FEsolving_execElem ( 2 ) ) > num % subStepMinHomog ) )
2009-05-07 21:57:36 +05:30
2019-06-15 20:12:16 +05:30
!$OMP PARALLEL DO PRIVATE(myNgrains)
elementLooping1 : do e = FEsolving_execElem ( 1 ) , FEsolving_execElem ( 2 )
myNgrains = homogenization_Ngrains ( material_homogenizationAt ( e ) )
2020-01-25 13:54:42 +05:30
IpLooping1 : do i = FEsolving_execIP ( 1 ) , FEsolving_execIP ( 2 )
2015-10-15 00:06:19 +05:30
2020-04-15 16:39:05 +05:30
if ( converged ( i , e ) ) then
2017-10-03 18:50:53 +05:30
#ifdef DEBUG
2020-07-02 04:55:24 +05:30
if ( debugHomog % extensive &
. and . ( ( e == debugHomog % element . and . i == debugHomog % ip ) &
. or . . not . debugHomog % selective ) ) then
2019-06-15 20:12:16 +05:30
write ( 6 , '(a,1x,f12.8,1x,a,1x,f12.8,1x,a,i8,1x,i2/)' ) '<< HOMOG >> winding forward from' , &
2020-04-15 16:39:05 +05:30
subFrac ( i , e ) , 'to current subFrac' , &
subFrac ( i , e ) + subStep ( i , e ) , 'in materialpoint_stressAndItsTangent at el ip' , e , i
2019-06-15 20:12:16 +05:30
endif
2011-03-29 12:57:19 +05:30
#endif
2015-10-15 00:06:19 +05:30
!---------------------------------------------------------------------------------------------------
2013-01-29 15:58:01 +05:30
! calculate new subStep and new subFrac
2020-04-15 16:39:05 +05:30
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
crystallite_partionedF0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = crystallite_partionedF ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
crystallite_partionedFp0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = crystallite_Fp ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
crystallite_partionedLp0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = crystallite_Lp ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
crystallite_partionedFi0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = crystallite_Fi ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
crystallite_partionedLi0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = crystallite_Li ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
crystallite_partionedS0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = crystallite_S ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
2019-06-15 20:12:16 +05:30
do g = 1 , myNgrains
plasticState ( material_phaseAt ( g , e ) ) % partionedState0 ( : , material_phasememberAt ( g , i , e ) ) = &
plasticState ( material_phaseAt ( g , e ) ) % state ( : , material_phasememberAt ( g , i , e ) )
do mySource = 1 , phase_Nsources ( material_phaseAt ( g , e ) )
sourceState ( material_phaseAt ( g , e ) ) % p ( mySource ) % partionedState0 ( : , material_phasememberAt ( g , i , e ) ) = &
sourceState ( material_phaseAt ( g , e ) ) % p ( mySource ) % state ( : , material_phasememberAt ( g , i , e ) )
enddo
enddo
if ( homogState ( material_homogenizationAt ( e ) ) % sizeState > 0 ) &
2020-01-23 17:18:18 +05:30
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 ( thermalState ( material_homogenizationAt ( e ) ) % sizeState > 0 ) &
2020-01-23 17:18:18 +05:30
thermalState ( material_homogenizationAt ( e ) ) % subState0 ( : , material_homogenizationMemberAt ( i , e ) ) = &
thermalState ( material_homogenizationAt ( e ) ) % State ( : , material_homogenizationMemberAt ( i , e ) )
2019-06-15 20:12:16 +05:30
if ( damageState ( material_homogenizationAt ( e ) ) % sizeState > 0 ) &
2020-01-23 17:18:18 +05:30
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
2015-05-28 22:32:23 +05:30
! cutback makes no sense
2019-06-15 20:12:16 +05:30
if ( . not . terminallyIll ) then ! so first signals terminally ill...
!$OMP CRITICAL (write2out)
2020-02-22 04:07:35 +05:30
write ( 6 , * ) 'Integration point ' , i , ' at element ' , e , ' terminally ill'
2019-06-15 20:12:16 +05:30
!$OMP END CRITICAL (write2out)
endif
2020-02-22 04:07:35 +05:30
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
2015-10-15 00:06:19 +05:30
2017-10-03 18:50:53 +05:30
#ifdef DEBUG
2020-07-02 04:55:24 +05:30
if ( debugHomog % extensive &
. and . ( ( e == debugHomog % element . and . i == debugHomog % ip ) &
. or . . not . debugHomog % selective ) ) then
2019-06-15 20:12:16 +05:30
write ( 6 , '(a,1x,f12.8,a,i8,1x,i2/)' ) &
2020-04-15 16:39:05 +05:30
'<< HOMOG >> cutback step in materialpoint_stressAndItsTangent with new subStep:' , &
subStep ( i , e ) , ' at el ip' , e , i
2019-06-15 20:12:16 +05:30
endif
2011-03-29 12:57:19 +05:30
#endif
2015-10-15 00:06:19 +05:30
2013-01-29 15:58:01 +05:30
!--------------------------------------------------------------------------------------------------
2020-04-15 12:23:25 +05:30
! restore
2020-04-15 16:39:05 +05:30
if ( subStep ( i , e ) < 1.0_pReal ) then ! protect against fake cutback from \Delta t = 2 to 1. Maybe that "trick" is not necessary anymore at all? I.e. start with \Delta t = 1
2020-04-15 12:23:25 +05:30
crystallite_Lp ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = crystallite_partionedLp0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
crystallite_Li ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = crystallite_partionedLi0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
2019-06-15 20:12:16 +05:30
endif ! maybe protecting everything from overwriting (not only L) makes even more sense
2020-04-15 12:23:25 +05:30
crystallite_Fp ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = crystallite_partionedFp0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
crystallite_Fi ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = crystallite_partionedFi0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
crystallite_S ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = crystallite_partionedS0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
2019-06-15 20:12:16 +05:30
do g = 1 , myNgrains
plasticState ( material_phaseAt ( g , e ) ) % state ( : , material_phasememberAt ( g , i , e ) ) = &
plasticState ( material_phaseAt ( g , e ) ) % partionedState0 ( : , material_phasememberAt ( g , i , e ) )
do mySource = 1 , phase_Nsources ( material_phaseAt ( g , e ) )
sourceState ( material_phaseAt ( g , e ) ) % p ( mySource ) % state ( : , material_phasememberAt ( g , i , e ) ) = &
sourceState ( material_phaseAt ( g , e ) ) % p ( mySource ) % partionedState0 ( : , material_phasememberAt ( g , i , e ) )
enddo
enddo
if ( homogState ( material_homogenizationAt ( e ) ) % sizeState > 0 ) &
2020-01-23 17:18:18 +05:30
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 ( thermalState ( material_homogenizationAt ( e ) ) % sizeState > 0 ) &
2020-01-23 17:18:18 +05:30
thermalState ( material_homogenizationAt ( e ) ) % State ( : , material_homogenizationMemberAt ( i , e ) ) = &
thermalState ( material_homogenizationAt ( e ) ) % subState0 ( : , material_homogenizationMemberAt ( i , e ) )
2019-06-15 20:12:16 +05:30
if ( damageState ( material_homogenizationAt ( e ) ) % sizeState > 0 ) &
2020-01-23 17:18:18 +05:30
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
2009-05-07 21:57:36 +05:30
2013-01-29 15:58:01 +05:30
!--------------------------------------------------------------------------------------------------
! deformation partitioning
2020-04-16 01:26:20 +05:30
!$OMP PARALLEL DO PRIVATE(myNgrains)
2019-06-15 20:12:16 +05:30
elementLooping2 : do e = FEsolving_execElem ( 1 ) , FEsolving_execElem ( 2 )
myNgrains = homogenization_Ngrains ( material_homogenizationAt ( e ) )
2020-01-25 13:54:42 +05:30
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
2020-04-16 01:26:20 +05:30
call partitionDeformation ( materialpoint_F0 ( 1 : 3 , 1 : 3 , i , e ) &
+ ( materialpoint_F ( 1 : 3 , 1 : 3 , i , e ) - materialpoint_F0 ( 1 : 3 , 1 : 3 , i , e ) ) &
* ( 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
2015-10-15 00:06:19 +05:30
2013-01-29 15:58:01 +05:30
!--------------------------------------------------------------------------------------------------
! crystallite integration
2020-04-15 16:39:05 +05:30
converged = crystallite_stress ( ) !ToDo: MD not sure if that is the best logic
2009-05-07 21:57:36 +05:30
2013-01-29 15:58:01 +05:30
!--------------------------------------------------------------------------------------------------
! state update
2020-04-16 01:26:20 +05:30
!$OMP PARALLEL DO
2019-06-15 20:12:16 +05:30
elementLooping3 : do e = FEsolving_execElem ( 1 ) , FEsolving_execElem ( 2 )
2020-01-25 13:54:42 +05:30
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
2020-04-16 01:26:20 +05:30
doneAndHappy ( 1 : 2 , i , e ) = updateState ( dt * subStep ( i , e ) , &
materialpoint_F0 ( 1 : 3 , 1 : 3 , i , e ) &
+ ( materialpoint_F ( 1 : 3 , 1 : 3 , i , e ) - materialpoint_F0 ( 1 : 3 , 1 : 3 , i , e ) ) &
* ( 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 ( updateJaco ) call crystallite_stressTangent
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 )
2020-01-25 13:54:42 +05:30
IpLooping4 : do i = FEsolving_execIP ( 1 ) , FEsolving_execIP ( 2 )
2019-06-15 20:12:16 +05:30
call averageStressAndItsTangent ( i , e )
enddo IpLooping4
enddo elementLooping4
!$OMP END PARALLEL DO
else
write ( 6 , '(/,a,/)' ) '<< HOMOG >> Material Point terminally ill'
endif
2015-10-15 00:06:19 +05:30
2012-08-10 21:28:17 +05:30
end subroutine materialpoint_stressAndItsTangent
2009-05-07 21:57:36 +05:30
2012-08-10 21:28:17 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief partition material point def grad onto constituents
!--------------------------------------------------------------------------------------------------
2020-04-14 11:15:39 +05:30
subroutine partitionDeformation ( subF , ip , el )
2009-05-07 21:57:36 +05:30
2020-04-14 13:19:03 +05:30
real ( pReal ) , intent ( in ) , dimension ( 3 , 3 ) :: &
subF
integer , intent ( in ) :: &
ip , & !< integration point
el !< element number
2020-07-03 20:15:11 +05:30
2020-04-14 13:19:03 +05:30
chosenHomogenization : select case ( homogenization_type ( material_homogenizationAt ( el ) ) )
2014-03-14 04:50:50 +05:30
2020-04-14 13:19:03 +05:30
case ( HOMOGENIZATION_NONE_ID ) chosenHomogenization
crystallite_partionedF ( 1 : 3 , 1 : 3 , 1 , ip , el ) = subF
2014-03-14 04:50:50 +05:30
2020-04-14 13:19:03 +05:30
case ( HOMOGENIZATION_ISOSTRAIN_ID ) chosenHomogenization
call mech_isostrain_partitionDeformation ( &
crystallite_partionedF ( 1 : 3 , 1 : 3 , 1 : homogenization_Ngrains ( material_homogenizationAt ( el ) ) , ip , el ) , &
subF )
2018-11-03 21:10:17 +05:30
2020-04-14 13:19:03 +05:30
case ( HOMOGENIZATION_RGC_ID ) chosenHomogenization
call mech_RGC_partitionDeformation ( &
crystallite_partionedF ( 1 : 3 , 1 : 3 , 1 : homogenization_Ngrains ( material_homogenizationAt ( el ) ) , ip , el ) , &
subF , &
ip , &
2020-07-02 01:50:22 +05:30
el )
2020-04-14 13:19:03 +05:30
end select chosenHomogenization
2009-05-07 21:57:36 +05:30
2019-01-13 14:18:47 +05:30
end subroutine partitionDeformation
2009-05-07 21:57:36 +05:30
2012-08-10 21:28:17 +05:30
!--------------------------------------------------------------------------------------------------
2015-10-15 00:06:19 +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 )
2013-11-27 13:34:05 +05:30
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
logical , dimension ( 2 ) :: updateState
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
updateState = &
updateState . and . &
mech_RGC_updateState ( crystallite_P ( 1 : 3 , 1 : 3 , 1 : homogenization_Ngrains ( material_homogenizationAt ( el ) ) , ip , el ) , &
crystallite_partionedF ( 1 : 3 , 1 : 3 , 1 : homogenization_Ngrains ( material_homogenizationAt ( el ) ) , ip , el ) , &
crystallite_partionedF0 ( 1 : 3 , 1 : 3 , 1 : homogenization_Ngrains ( material_homogenizationAt ( el ) ) , ip , el ) , &
subF , &
subdt , &
crystallite_dPdF ( 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 , 1 : homogenization_Ngrains ( material_homogenizationAt ( el ) ) , ip , el ) , &
ip , &
2020-07-02 01:50:22 +05:30
el )
2020-04-14 13:19:03 +05:30
end select chosenHomogenization
chosenThermal : select case ( thermal_type ( material_homogenizationAt ( el ) ) )
case ( THERMAL_adiabatic_ID ) chosenThermal
updateState = &
updateState . and . &
thermal_adiabatic_updateState ( subdt , &
ip , &
el )
end select chosenThermal
chosenDamage : select case ( damage_type ( material_homogenizationAt ( el ) ) )
case ( DAMAGE_local_ID ) chosenDamage
updateState = &
updateState . and . &
damage_local_updateState ( subdt , &
ip , &
el )
end select chosenDamage
2015-05-28 22:32:23 +05:30
2019-01-13 14:18:47 +05:30
end function updateState
2009-05-07 21:57:36 +05:30
2012-08-10 21:28:17 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief derive average stress and stiffness from constituent quantities
!--------------------------------------------------------------------------------------------------
2019-01-13 14:18:47 +05:30
subroutine averageStressAndItsTangent ( ip , el )
2009-07-22 21:37:19 +05:30
2020-04-14 13:19:03 +05:30
integer , intent ( in ) :: &
ip , & !< integration point
el !< element number
chosenHomogenization : select case ( homogenization_type ( material_homogenizationAt ( el ) ) )
case ( HOMOGENIZATION_NONE_ID ) chosenHomogenization
materialpoint_P ( 1 : 3 , 1 : 3 , ip , el ) = crystallite_P ( 1 : 3 , 1 : 3 , 1 , ip , el )
materialpoint_dPdF ( 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 , ip , el ) = crystallite_dPdF ( 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 , 1 , ip , el )
case ( HOMOGENIZATION_ISOSTRAIN_ID ) chosenHomogenization
call mech_isostrain_averageStressAndItsTangent ( &
materialpoint_P ( 1 : 3 , 1 : 3 , ip , el ) , &
materialpoint_dPdF ( 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 , ip , el ) , &
crystallite_P ( 1 : 3 , 1 : 3 , 1 : homogenization_Ngrains ( material_homogenizationAt ( el ) ) , ip , el ) , &
crystallite_dPdF ( 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 , 1 : homogenization_Ngrains ( material_homogenizationAt ( el ) ) , ip , el ) , &
homogenization_typeInstance ( material_homogenizationAt ( el ) ) )
case ( HOMOGENIZATION_RGC_ID ) chosenHomogenization
call mech_RGC_averageStressAndItsTangent ( &
materialpoint_P ( 1 : 3 , 1 : 3 , ip , el ) , &
materialpoint_dPdF ( 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 , ip , el ) , &
crystallite_P ( 1 : 3 , 1 : 3 , 1 : homogenization_Ngrains ( material_homogenizationAt ( el ) ) , ip , el ) , &
crystallite_dPdF ( 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 , 1 : homogenization_Ngrains ( material_homogenizationAt ( el ) ) , ip , el ) , &
homogenization_typeInstance ( material_homogenizationAt ( el ) ) )
end select chosenHomogenization
2009-07-22 21:37:19 +05:30
2019-01-13 14:18:47 +05:30
end subroutine averageStressAndItsTangent
2009-07-22 21:37:19 +05:30
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
2019-12-10 11:44:39 +05:30
character ( len = pStringLen ) :: group_base , group
2020-03-17 05:09:32 +05:30
2019-09-22 19:23:03 +05:30
!real(pReal), dimension(:,:,:), allocatable :: temp
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/materialpoint/' / / 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
2019-12-10 11:44:39 +05:30
group = trim ( group_base ) / / '/generic'
2019-12-19 12:19:53 +05:30
call results_closeGroup ( results_addGroup ( group ) )
2019-07-16 05:38:18 +05:30
!temp = reshape(materialpoint_F,[3,3,discretization_nIP*discretization_nElem])
!call results_writeDataset(group,temp,'F',&
2020-03-17 05:09:32 +05:30
! 'deformation gradient','1')
2019-07-16 05:38:18 +05:30
!temp = reshape(materialpoint_P,[3,3,discretization_nIP*discretization_nElem])
!call results_writeDataset(group,temp,'P',&
2020-08-13 00:44:00 +05:30
! '1st Piola-Kirchhoff stress','Pa')
2020-03-17 05:09:32 +05:30
2019-12-10 11:44:39 +05:30
group = trim ( group_base ) / / '/mech'
call results_closeGroup ( results_addGroup ( group ) )
2019-04-30 22:15:16 +05:30
select case ( material_homogenization_type ( p ) )
2019-05-01 02:35:21 +05:30
case ( HOMOGENIZATION_rgc_ID )
call mech_RGC_results ( homogenization_typeInstance ( p ) , group )
2019-04-30 22:15:16 +05:30
end select
2019-07-16 05:38:18 +05:30
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_LOCAL_ID )
call damage_local_results ( p , group )
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_ADIABATIC_ID )
call thermal_adiabatic_results ( p , group )
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
2012-08-10 21:28:17 +05:30
end module homogenization