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-05-15 02:42:32 +05:30
use prec
2019-05-18 10:24:45 +05:30
use IO
use config
use debug
use math
2019-05-01 02:23:32 +05:30
use material
2019-05-18 10:24:45 +05:30
use numerics
use constitutive
use crystallite
use FEsolving
2019-06-07 02:19:17 +05:30
use mesh
use discretization
2019-05-28 15:36:21 +05:30
use thermal_isothermal
use thermal_adiabatic
use thermal_conduction
use damage_none
use damage_local
use damage_nonlocal
2019-05-18 10:24:45 +05:30
use results
use HDF5_utilities
2014-08-21 23:18:20 +05:30
2009-05-07 21:57:36 +05:30
implicit none
2013-01-29 15:58:01 +05:30
private
2019-06-11 13:18:07 +05:30
!--------------------------------------------------------------------------------------------------
! General variables for the homogenization at a material point
real ( pReal ) , dimension ( : , : , : , : ) , allocatable , public :: &
2013-01-29 15:58:01 +05:30
materialpoint_F0 , & !< def grad of IP at start of FE increment
materialpoint_F , & !< def grad of IP to be reached at end of FE increment
materialpoint_P !< first P--K stress of IP
real ( pReal ) , dimension ( : , : , : , : , : , : ) , allocatable , public :: &
materialpoint_dPdF !< tangent of first P--K stress at IP
real ( pReal ) , dimension ( : , : , : ) , allocatable , public :: &
materialpoint_results !< results array of material point
2019-04-06 00:15:56 +05:30
integer , public , protected :: &
2013-01-29 15:58:01 +05:30
materialpoint_sizeResults , &
2015-05-28 22:32:23 +05:30
thermal_maxSizePostResults , &
2018-12-22 13:30:57 +05:30
damage_maxSizePostResults
2013-01-29 15:58:01 +05:30
2019-06-11 13:18:07 +05:30
real ( pReal ) , dimension ( : , : , : , : ) , allocatable :: &
2013-01-29 15:58:01 +05:30
materialpoint_subF0 , & !< def grad of IP at beginning of homogenization increment
materialpoint_subF !< def grad of IP to be reached at end of homog inc
2019-06-11 13:18:07 +05:30
real ( pReal ) , dimension ( : , : ) , allocatable :: &
2013-01-29 15:58:01 +05:30
materialpoint_subFrac , &
materialpoint_subStep , &
2013-10-19 02:26:10 +05:30
materialpoint_subdt
2019-06-11 13:18:07 +05:30
logical , dimension ( : , : ) , allocatable :: &
2013-01-29 15:58:01 +05:30
materialpoint_requested , &
materialpoint_converged
2019-06-11 13:18:07 +05:30
logical , dimension ( : , : , : ) , allocatable :: &
2013-01-29 15:58:01 +05:30
materialpoint_doneAndHappy
2019-04-06 01:12:07 +05:30
interface
module subroutine mech_none_init
end subroutine mech_none_init
module subroutine mech_isostrain_init
end subroutine mech_isostrain_init
2019-05-18 10:53:46 +05:30
module subroutine mech_RGC_init
end subroutine mech_RGC_init
2019-04-06 01:12:07 +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
2019-05-18 10:53:46 +05:30
module subroutine mech_RGC_partitionDeformation ( F , avgF , instance , of )
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
2019-04-06 01:12:07 +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
real ( pReal ) , dimension ( : , : , : ) , intent ( in ) :: P !< partitioned stresses
real ( pReal ) , dimension ( : , : , : , : , : ) , intent ( in ) :: dPdF !< partitioned stiffnesses
integer , intent ( in ) :: instance
end subroutine mech_isostrain_averageStressAndItsTangent
2019-05-18 10:53:46 +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
real ( pReal ) , dimension ( : , : , : ) , intent ( in ) :: P !< partitioned stresses
real ( pReal ) , dimension ( : , : , : , : , : ) , intent ( in ) :: dPdF !< partitioned stiffnesses
integer , intent ( in ) :: instance
end subroutine mech_RGC_averageStressAndItsTangent
module function mech_RGC_updateState ( P , F , F0 , avgF , dt , dPdF , ip , el )
2019-05-18 11:09:55 +05:30
logical , dimension ( 2 ) :: mech_RGC_updateState
real ( pReal ) , dimension ( : , : , : ) , intent ( in ) :: &
P , & !< partitioned stresses
F , & !< partitioned deformation gradients
F0 !< partitioned initial deformation gradients
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 ) :: &
ip , & !< integration point number
el !< element number
end function mech_RGC_updateState
2019-05-18 10:53:46 +05:30
2019-04-06 01:12:07 +05:30
2019-05-18 11:09:55 +05:30
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
2019-04-06 01:12:07 +05:30
end interface
2013-01-29 15:58:01 +05:30
public :: &
homogenization_init , &
materialpoint_stressAndItsTangent , &
2019-04-30 22:15:16 +05:30
materialpoint_postResults , &
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
2019-04-06 00:15:56 +05:30
integer , parameter :: FILEUNIT = 200
integer :: e , i , p
integer , dimension ( : , : ) , pointer :: thisSize
integer , dimension ( : ) , pointer :: thisNoutput
2012-08-10 21:28:17 +05:30
character ( len = 64 ) , dimension ( : , : ) , pointer :: thisOutput
2013-11-27 13:34:05 +05:30
character ( len = 32 ) :: outputName !< name of output, intermediate fix until HDF5 output is ready
2018-02-16 20:06:18 +05:30
logical :: valid
2014-09-26 16:04:36 +05:30
2019-04-06 01:12:07 +05:30
if ( any ( homogenization_type == HOMOGENIZATION_NONE_ID ) ) call mech_none_init
if ( any ( homogenization_type == HOMOGENIZATION_ISOSTRAIN_ID ) ) call mech_isostrain_init
2019-05-18 10:53:46 +05:30
if ( any ( homogenization_type == HOMOGENIZATION_RGC_ID ) ) call mech_RGC_init
2015-09-24 23:20:11 +05:30
2019-02-13 03:22:21 +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-03-09 03:46:56 +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
2012-08-10 21:28:17 +05:30
!--------------------------------------------------------------------------------------------------
! write description file for homogenization output
2019-05-17 10:06:30 +05:30
mainProcess : if ( worldrank == 0 ) then
2015-03-25 21:32:30 +05:30
call IO_write_jobFile ( FILEUNIT , 'outputHomogenization' )
2018-06-27 00:24:54 +05:30
do p = 1 , size ( config_homogenization )
2019-03-10 15:32:32 +05:30
if ( any ( material_homogenizationAt == p ) ) then
2019-06-15 17:56:32 +05:30
write ( FILEUNIT , '(/,a,/)' ) '[' / / trim ( config_name_homogenization ( p ) ) / / ']'
2019-05-17 10:06:30 +05:30
write ( FILEUNIT , '(a)' ) '(type) n/a'
write ( FILEUNIT , '(a,i4)' ) '(ngrains)' / / char ( 9 ) , homogenization_Ngrains ( p )
2015-05-28 22:32:23 +05:30
i = thermal_typeInstance ( p ) ! which instance of this thermal type
2018-02-16 20:06:18 +05:30
valid = . true . ! assume valid
2015-05-28 22:32:23 +05:30
select case ( thermal_type ( p ) ) ! split per thermal type
case ( THERMAL_isothermal_ID )
outputName = THERMAL_isothermal_label
thisNoutput = > null ( )
thisOutput = > null ( )
thisSize = > null ( )
case ( THERMAL_adiabatic_ID )
outputName = THERMAL_adiabatic_label
thisNoutput = > thermal_adiabatic_Noutput
thisOutput = > thermal_adiabatic_output
thisSize = > thermal_adiabatic_sizePostResult
case ( THERMAL_conduction_ID )
outputName = THERMAL_conduction_label
thisNoutput = > thermal_conduction_Noutput
thisOutput = > thermal_conduction_output
thisSize = > thermal_conduction_sizePostResult
case default
2018-02-16 20:06:18 +05:30
valid = . false .
2015-10-15 00:06:19 +05:30
end select
2018-02-16 20:06:18 +05:30
if ( valid ) then
2015-05-28 22:32:23 +05:30
write ( FILEUNIT , '(a)' ) '(thermal)' / / char ( 9 ) / / trim ( outputName )
if ( thermal_type ( p ) / = THERMAL_isothermal_ID ) then
do e = 1 , thisNoutput ( i )
write ( FILEUNIT , '(a,i4)' ) trim ( thisOutput ( e , i ) ) / / char ( 9 ) , thisSize ( e , i )
enddo
2015-10-15 00:06:19 +05:30
endif
endif
2019-05-17 10:06:30 +05:30
2015-05-28 22:32:23 +05:30
i = damage_typeInstance ( p ) ! which instance of this damage type
2018-02-16 20:06:18 +05:30
valid = . true . ! assume valid
2015-05-28 22:32:23 +05:30
select case ( damage_type ( p ) ) ! split per damage type
case ( DAMAGE_none_ID )
outputName = DAMAGE_none_label
thisNoutput = > null ( )
thisOutput = > null ( )
thisSize = > null ( )
case ( DAMAGE_local_ID )
outputName = DAMAGE_local_label
thisNoutput = > damage_local_Noutput
thisOutput = > damage_local_output
thisSize = > damage_local_sizePostResult
case ( DAMAGE_nonlocal_ID )
outputName = DAMAGE_nonlocal_label
thisNoutput = > damage_nonlocal_Noutput
thisOutput = > damage_nonlocal_output
thisSize = > damage_nonlocal_sizePostResult
case default
2018-02-16 20:06:18 +05:30
valid = . false .
2015-10-15 00:06:19 +05:30
end select
2018-02-16 20:06:18 +05:30
if ( valid ) then
2015-05-28 22:32:23 +05:30
write ( FILEUNIT , '(a)' ) '(damage)' / / char ( 9 ) / / trim ( outputName )
if ( damage_type ( p ) / = DAMAGE_none_ID ) then
do e = 1 , thisNoutput ( i )
write ( FILEUNIT , '(a,i4)' ) trim ( thisOutput ( e , i ) ) / / char ( 9 ) , thisSize ( e , i )
enddo
2015-10-15 00:06:19 +05:30
endif
endif
2015-04-21 19:40:34 +05:30
endif
2014-09-26 16:04:36 +05:30
enddo
2015-03-25 21:32:30 +05:30
close ( FILEUNIT )
2019-05-17 10:06:30 +05:30
endif mainProcess
2014-09-22 23:45:19 +05:30
2018-06-27 00:20:06 +05:30
call config_deallocate ( 'material.config/homogenization' )
2012-08-10 21:28:17 +05:30
!--------------------------------------------------------------------------------------------------
! allocate and initialize global variables
2019-06-07 02:19:17 +05:30
allocate ( materialpoint_dPdF ( 3 , 3 , 3 , 3 , discretization_nIP , discretization_nElem ) , source = 0.0_pReal )
allocate ( materialpoint_F0 ( 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
allocate ( materialpoint_F ( 3 , 3 , discretization_nIP , discretization_nElem ) , source = 0.0_pReal )
2013-12-16 17:28:03 +05:30
materialpoint_F = materialpoint_F0 ! initialize to identity
2019-06-07 02:19:17 +05:30
allocate ( materialpoint_subF0 ( 3 , 3 , discretization_nIP , discretization_nElem ) , source = 0.0_pReal )
allocate ( materialpoint_subF ( 3 , 3 , discretization_nIP , discretization_nElem ) , source = 0.0_pReal )
allocate ( materialpoint_P ( 3 , 3 , discretization_nIP , discretization_nElem ) , source = 0.0_pReal )
allocate ( materialpoint_subFrac ( discretization_nIP , discretization_nElem ) , source = 0.0_pReal )
allocate ( materialpoint_subStep ( discretization_nIP , discretization_nElem ) , source = 0.0_pReal )
allocate ( materialpoint_subdt ( discretization_nIP , discretization_nElem ) , source = 0.0_pReal )
allocate ( materialpoint_requested ( discretization_nIP , discretization_nElem ) , source = . false . )
allocate ( materialpoint_converged ( discretization_nIP , discretization_nElem ) , source = . true . )
allocate ( materialpoint_doneAndHappy ( 2 , discretization_nIP , discretization_nElem ) , source = . true . )
2013-12-16 17:28:03 +05:30
2012-08-10 21:28:17 +05:30
!--------------------------------------------------------------------------------------------------
2013-10-16 18:34:59 +05:30
! allocate and initialize global state and postresutls variables
2019-04-06 00:15:56 +05:30
thermal_maxSizePostResults = 0
damage_maxSizePostResults = 0
2018-06-27 00:24:54 +05:30
do p = 1 , size ( config_homogenization )
2015-05-28 22:32:23 +05:30
thermal_maxSizePostResults = max ( thermal_maxSizePostResults , thermalState ( p ) % sizePostResults )
damage_maxSizePostResults = max ( damage_maxSizePostResults , damageState ( p ) % sizePostResults )
2015-10-15 00:06:19 +05:30
enddo
2015-03-30 15:15:10 +05:30
2013-01-29 15:58:01 +05:30
materialpoint_sizeResults = 1 & ! grain count
2019-05-18 11:09:55 +05:30
+ 1 + thermal_maxSizePostResults &
2015-05-28 22:32:23 +05:30
+ damage_maxSizePostResults &
2013-01-29 15:58:01 +05:30
+ homogenization_maxNgrains * ( 1 + crystallite_maxSizePostResults & ! crystallite size & crystallite results
2015-05-28 22:32:23 +05:30
+ 1 + constitutive_plasticity_maxSizePostResults & ! constitutive size & constitutive results
2015-10-15 00:06:19 +05:30
+ constitutive_source_maxSizePostResults )
2019-06-07 02:19:17 +05:30
allocate ( materialpoint_results ( materialpoint_sizeResults , discretization_nIP , discretization_nElem ) )
2015-10-15 00:06:19 +05:30
2018-06-02 22:58:08 +05:30
write ( 6 , '(/,a)' ) ' <<<+- homogenization init -+>>>'
2014-10-10 01:53:06 +05:30
2019-04-06 00:15:56 +05:30
if ( iand ( debug_level ( debug_homogenization ) , debug_levelBasic ) / = 0 ) then
2013-10-16 18:34:59 +05:30
write ( 6 , '(a32,1x,7(i8,1x))' ) 'materialpoint_dPdF: ' , shape ( materialpoint_dPdF )
write ( 6 , '(a32,1x,7(i8,1x))' ) 'materialpoint_F0: ' , shape ( materialpoint_F0 )
write ( 6 , '(a32,1x,7(i8,1x))' ) 'materialpoint_F: ' , shape ( materialpoint_F )
write ( 6 , '(a32,1x,7(i8,1x))' ) 'materialpoint_subF0: ' , shape ( materialpoint_subF0 )
write ( 6 , '(a32,1x,7(i8,1x))' ) 'materialpoint_subF: ' , shape ( materialpoint_subF )
write ( 6 , '(a32,1x,7(i8,1x))' ) 'materialpoint_P: ' , shape ( materialpoint_P )
write ( 6 , '(a32,1x,7(i8,1x))' ) 'materialpoint_subFrac: ' , shape ( materialpoint_subFrac )
write ( 6 , '(a32,1x,7(i8,1x))' ) 'materialpoint_subStep: ' , shape ( materialpoint_subStep )
write ( 6 , '(a32,1x,7(i8,1x))' ) 'materialpoint_subdt: ' , shape ( materialpoint_subdt )
write ( 6 , '(a32,1x,7(i8,1x))' ) 'materialpoint_requested: ' , shape ( materialpoint_requested )
write ( 6 , '(a32,1x,7(i8,1x))' ) 'materialpoint_converged: ' , shape ( materialpoint_converged )
write ( 6 , '(a32,1x,7(i8,1x),/)' ) 'materialpoint_doneAndHappy: ' , shape ( materialpoint_doneAndHappy )
2013-01-29 15:58:01 +05:30
endif
flush ( 6 )
2015-10-15 00:06:19 +05:30
2019-06-07 02:19:17 +05:30
if ( debug_g < 1 . or . debug_g > homogenization_Ngrains ( material_homogenizationAt ( debug_e ) ) ) &
2019-04-06 00:15:56 +05:30
call IO_error ( 602 , ext_msg = 'constituent' , el = debug_e , g = debug_g )
2015-10-15 00:06:19 +05:30
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
2012-08-10 21:28:17 +05:30
real ( pReal ) , intent ( in ) :: dt !< time increment
logical , intent ( in ) :: updateJaco !< initiating Jacobian update
2019-04-06 00:15:56 +05:30
integer :: &
2013-01-29 15:58:01 +05:30
NiterationHomog , &
NiterationMPstate , &
g , & !< grain number
i , & !< integration point number
e , & !< element number
2015-05-28 22:32:23 +05:30
mySource , &
2013-01-29 15:58:01 +05:30
myNgrains
2009-05-07 21:57:36 +05:30
2018-09-20 10:57:12 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
if ( iand ( debug_level ( debug_homogenization ) , debug_levelBasic ) / = 0 ) then
2013-10-16 18:34:59 +05:30
write ( 6 , '(/a,i5,1x,i2)' ) '<< HOMOG >> Material Point start at el ip ' , debug_e , debug_i
2012-10-02 18:23:25 +05:30
write ( 6 , '(a,/,3(12x,3(f14.9,1x)/))' ) '<< HOMOG >> F0' , &
2018-09-20 01:15:12 +05:30
transpose ( materialpoint_F0 ( 1 : 3 , 1 : 3 , debug_i , debug_e ) )
2012-10-02 18:23:25 +05:30
write ( 6 , '(a,/,3(12x,3(f14.9,1x)/))' ) '<< HOMOG >> F' , &
2018-09-20 01:15:12 +05:30
transpose ( materialpoint_F ( 1 : 3 , 1 : 3 , debug_i , debug_e ) )
2010-03-19 19:44:08 +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
!--------------------------------------------------------------------------------------------------
2012-11-07 21:13:29 +05:30
! initialize restoration points of ...
do e = FEsolving_execElem ( 1 ) , FEsolving_execElem ( 2 )
2019-06-07 00:44:37 +05:30
myNgrains = homogenization_Ngrains ( material_homogenizationAt ( e ) )
2019-05-16 21:54:54 +05:30
do i = FEsolving_execIP ( 1 , e ) , FEsolving_execIP ( 2 , e ) ;
do g = 1 , myNgrains
2019-06-14 12:39:16 +05:30
plasticState ( material_phaseAt ( g , e ) ) % partionedState0 ( : , material_phasememberAt ( g , i , e ) ) = &
plasticState ( material_phaseAt ( g , e ) ) % state0 ( : , material_phasememberAt ( g , i , e ) )
2019-06-14 12:32:28 +05:30
do mySource = 1 , phase_Nsources ( material_phaseAt ( g , e ) )
2019-06-14 12:39:16 +05:30
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 ) )
2019-05-16 21:54:54 +05:30
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 )
2014-09-03 01:46:33 +05:30
2015-10-15 00:06:19 +05:30
enddo
2014-09-03 01:46:33 +05:30
2014-05-22 17:37:50 +05:30
2019-05-16 21:54:54 +05:30
materialpoint_subF0 ( 1 : 3 , 1 : 3 , i , e ) = materialpoint_F0 ( 1 : 3 , 1 : 3 , i , e )
2009-05-07 21:57:36 +05:30
materialpoint_subFrac ( i , e ) = 0.0_pReal
2009-11-10 19:06:27 +05:30
materialpoint_subStep ( i , e ) = 1.0_pReal / subStepSizeHomog ! <<added to adopt flexibility in cutback size>>
2009-06-16 14:33:30 +05:30
materialpoint_converged ( i , e ) = . false . ! pretend failed step of twice the required size
materialpoint_requested ( i , e ) = . true . ! everybody requires calculation
2019-05-16 21:54:54 +05:30
if ( homogState ( material_homogenizationAt ( e ) ) % sizeState > 0 ) &
homogState ( material_homogenizationAt ( e ) ) % subState0 ( : , mappingHomogenization ( 1 , i , e ) ) = &
homogState ( material_homogenizationAt ( e ) ) % State0 ( : , mappingHomogenization ( 1 , i , e ) ) ! ...internal homogenization state
if ( thermalState ( material_homogenizationAt ( e ) ) % sizeState > 0 ) &
thermalState ( material_homogenizationAt ( e ) ) % subState0 ( : , mappingHomogenization ( 1 , i , e ) ) = &
thermalState ( material_homogenizationAt ( e ) ) % State0 ( : , mappingHomogenization ( 1 , i , e ) ) ! ...internal thermal state
if ( damageState ( material_homogenizationAt ( e ) ) % sizeState > 0 ) &
damageState ( material_homogenizationAt ( e ) ) % subState0 ( : , mappingHomogenization ( 1 , i , e ) ) = &
damageState ( material_homogenizationAt ( e ) ) % State0 ( : , mappingHomogenization ( 1 , i , e ) ) ! ...internal damage state
enddo
2009-05-07 21:57:36 +05:30
enddo
2019-05-16 21:54:54 +05:30
2019-04-06 00:15:56 +05:30
NiterationHomog = 0
2015-10-15 00:06:19 +05:30
2013-01-29 15:58:01 +05:30
cutBackLooping : do while ( . not . terminallyIll . and . &
any ( materialpoint_subStep ( : , FEsolving_execELem ( 1 ) : FEsolving_execElem ( 2 ) ) > subStepMinHomog ) )
2009-05-07 21:57:36 +05:30
2010-11-03 22:52:48 +05:30
!$OMP PARALLEL DO PRIVATE(myNgrains)
2013-01-29 15:58:01 +05:30
elementLooping1 : do e = FEsolving_execElem ( 1 ) , FEsolving_execElem ( 2 )
2019-06-07 00:44:37 +05:30
myNgrains = homogenization_Ngrains ( material_homogenizationAt ( e ) )
2013-01-29 15:58:01 +05:30
IpLooping1 : do i = FEsolving_execIP ( 1 , e ) , FEsolving_execIP ( 2 , e )
2015-10-15 00:06:19 +05:30
2019-05-16 21:54:54 +05:30
converged : if ( materialpoint_converged ( i , e ) ) then
2017-10-03 18:50:53 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
if ( iand ( debug_level ( debug_homogenization ) , debug_levelExtensive ) / = 0 &
2015-10-15 00:06:19 +05:30
. and . ( ( e == debug_e . and . i == debug_i ) &
2019-04-06 00:15:56 +05:30
. or . . not . iand ( debug_level ( debug_homogenization ) , debug_levelSelective ) / = 0 ) ) then
2012-11-21 22:28:14 +05:30
write ( 6 , '(a,1x,f12.8,1x,a,1x,f12.8,1x,a,i8,1x,i2/)' ) '<< HOMOG >> winding forward from' , &
2011-03-29 12:57:19 +05:30
materialpoint_subFrac ( i , e ) , 'to current materialpoint_subFrac' , &
2012-11-21 22:28:14 +05:30
materialpoint_subFrac ( i , e ) + materialpoint_subStep ( i , e ) , 'in materialpoint_stressAndItsTangent at el ip' , e , i
2009-08-24 13:46:01 +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
2009-08-27 17:40:06 +05:30
materialpoint_subFrac ( i , e ) = materialpoint_subFrac ( i , e ) + materialpoint_subStep ( i , e )
2009-11-10 19:06:27 +05:30
materialpoint_subStep ( i , e ) = min ( 1.0_pReal - materialpoint_subFrac ( i , e ) , &
2015-05-28 22:32:23 +05:30
stepIncreaseHomog * materialpoint_subStep ( i , e ) ) ! introduce flexibility for step increase/acceleration
2015-10-15 00:06:19 +05:30
2013-01-29 15:58:01 +05:30
steppingNeeded : if ( materialpoint_subStep ( i , e ) > subStepMinHomog ) then
2015-10-15 00:06:19 +05:30
2009-06-16 14:33:30 +05:30
! wind forward grain starting point of...
2019-05-14 23:22:48 +05:30
crystallite_partionedF0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = &
2019-05-16 21:54:54 +05:30
crystallite_partionedF ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
2015-10-15 00:06:19 +05:30
2019-05-14 23:22:48 +05:30
crystallite_partionedFp0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = &
2019-05-16 21:54:54 +05:30
crystallite_Fp ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
2015-10-15 00:06:19 +05:30
2019-05-14 23:22:48 +05:30
crystallite_partionedLp0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = &
2019-05-16 21:54:54 +05:30
crystallite_Lp ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
2015-10-15 00:06:19 +05:30
2019-05-14 23:22:48 +05:30
crystallite_partionedFi0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = &
2019-05-16 21:54:54 +05:30
crystallite_Fi ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
2015-10-15 00:06:19 +05:30
2019-05-14 23:22:48 +05:30
crystallite_partionedLi0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = &
2019-05-16 21:54:54 +05:30
crystallite_Li ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
2015-10-15 00:06:19 +05:30
2019-05-14 23:22:48 +05:30
crystallite_partionedS0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = &
2019-05-16 21:54:54 +05:30
crystallite_S ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
2015-10-15 00:06:19 +05:30
2015-05-28 22:32:23 +05:30
do g = 1 , myNgrains
2019-06-14 12:39:16 +05:30
plasticState ( material_phaseAt ( g , e ) ) % partionedState0 ( : , material_phasememberAt ( g , i , e ) ) = &
plasticState ( material_phaseAt ( g , e ) ) % state ( : , material_phasememberAt ( g , i , e ) )
2019-06-14 12:32:28 +05:30
do mySource = 1 , phase_Nsources ( material_phaseAt ( g , e ) )
2019-06-14 12:39:16 +05:30
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 ) )
2015-05-28 22:32:23 +05:30
enddo
2015-10-15 00:06:19 +05:30
enddo
2019-05-16 21:54:54 +05:30
if ( homogState ( material_homogenizationAt ( e ) ) % sizeState > 0 ) &
2019-03-10 15:53:39 +05:30
homogState ( material_homogenizationAt ( e ) ) % subState0 ( : , mappingHomogenization ( 1 , i , e ) ) = &
2019-05-16 21:54:54 +05:30
homogState ( material_homogenizationAt ( e ) ) % State ( : , mappingHomogenization ( 1 , i , e ) )
if ( thermalState ( material_homogenizationAt ( e ) ) % sizeState > 0 ) &
2019-03-10 15:53:39 +05:30
thermalState ( material_homogenizationAt ( e ) ) % subState0 ( : , mappingHomogenization ( 1 , i , e ) ) = &
2019-05-16 21:54:54 +05:30
thermalState ( material_homogenizationAt ( e ) ) % State ( : , mappingHomogenization ( 1 , i , e ) )
if ( damageState ( material_homogenizationAt ( e ) ) % sizeState > 0 ) &
2019-03-10 15:53:39 +05:30
damageState ( material_homogenizationAt ( e ) ) % subState0 ( : , mappingHomogenization ( 1 , i , e ) ) = &
2019-05-16 21:54:54 +05:30
damageState ( material_homogenizationAt ( e ) ) % State ( : , mappingHomogenization ( 1 , i , e ) )
materialpoint_subF0 ( 1 : 3 , 1 : 3 , i , e ) = materialpoint_subF ( 1 : 3 , 1 : 3 , i , e )
2013-01-29 15:58:01 +05:30
endif steppingNeeded
else converged
2019-05-16 21:54:54 +05:30
if ( ( myNgrains == 1 . and . materialpoint_subStep ( i , e ) < = 1.0 ) . or . & ! single grain already tried internal subStepping in crystallite
2015-05-28 22:32:23 +05:30
subStepSizeHomog * materialpoint_subStep ( i , e ) < = subStepMinHomog ) then ! would require too small subStep
! cutback makes no sense
2012-12-16 16:24:13 +05:30
!$OMP FLUSH(terminallyIll)
2015-05-28 22:32:23 +05:30
if ( . not . terminallyIll ) then ! so first signals terminally ill...
2012-10-18 19:18:06 +05:30
!$OMP CRITICAL (write2out)
write ( 6 , * ) 'Integration point ' , i , ' at element ' , e , ' terminally ill'
!$OMP END CRITICAL (write2out)
endif
2010-11-03 22:52:48 +05:30
!$OMP CRITICAL (setTerminallyIll)
2015-05-28 22:32:23 +05:30
terminallyIll = . true . ! ...and kills all others
2010-11-03 22:52:48 +05:30
!$OMP END CRITICAL (setTerminallyIll)
2015-05-28 22:32:23 +05:30
else ! cutback makes sense
materialpoint_subStep ( i , e ) = subStepSizeHomog * materialpoint_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
2019-04-06 00:15:56 +05:30
if ( iand ( debug_level ( debug_homogenization ) , debug_levelExtensive ) / = 0 &
2012-03-09 01:55:28 +05:30
. and . ( ( e == debug_e . and . i == debug_i ) &
2019-04-06 00:15:56 +05:30
. or . . not . iand ( debug_level ( debug_homogenization ) , debug_levelSelective ) / = 0 ) ) then
2012-11-21 22:28:14 +05:30
write ( 6 , '(a,1x,f12.8,a,i8,1x,i2/)' ) &
2011-08-02 18:06:08 +05:30
'<< HOMOG >> cutback step in materialpoint_stressAndItsTangent with new materialpoint_subStep:' , &
2015-10-15 00:06:19 +05:30
materialpoint_subStep ( i , e ) , ' at el ip' , e , i
2010-09-02 02:34:02 +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
!--------------------------------------------------------------------------------------------------
! restore...
2019-05-14 23:16:25 +05:30
if ( materialpoint_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
crystallite_Lp ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = &
2019-05-16 21:54:54 +05:30
crystallite_partionedLp0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
2019-05-14 23:16:25 +05:30
crystallite_Li ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = &
2019-05-16 21:54:54 +05:30
crystallite_partionedLi0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
2019-05-14 23:16:25 +05:30
endif ! maybe protecting everything from overwriting (not only L) makes even more sense
2015-05-28 22:32:23 +05:30
crystallite_Fp ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = &
2019-05-16 21:54:54 +05:30
crystallite_partionedFp0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
2015-05-28 22:32:23 +05:30
crystallite_Fi ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = &
2019-05-16 21:54:54 +05:30
crystallite_partionedFi0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
2019-03-09 05:03:11 +05:30
crystallite_S ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e ) = &
2019-05-16 21:54:54 +05:30
crystallite_partionedS0 ( 1 : 3 , 1 : 3 , 1 : myNgrains , i , e )
2015-05-28 22:32:23 +05:30
do g = 1 , myNgrains
2019-06-14 12:39:16 +05:30
plasticState ( material_phaseAt ( g , e ) ) % state ( : , material_phasememberAt ( g , i , e ) ) = &
plasticState ( material_phaseAt ( g , e ) ) % partionedState0 ( : , material_phasememberAt ( g , i , e ) )
2019-06-14 12:32:28 +05:30
do mySource = 1 , phase_Nsources ( material_phaseAt ( g , e ) )
2019-06-14 12:39:16 +05:30
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 ) )
2015-05-28 22:32:23 +05:30
enddo
2015-10-15 00:06:19 +05:30
enddo
2019-05-16 21:54:54 +05:30
if ( homogState ( material_homogenizationAt ( e ) ) % sizeState > 0 ) &
2019-03-10 15:53:39 +05:30
homogState ( material_homogenizationAt ( e ) ) % State ( : , mappingHomogenization ( 1 , i , e ) ) = &
2019-05-16 21:54:54 +05:30
homogState ( material_homogenizationAt ( e ) ) % subState0 ( : , mappingHomogenization ( 1 , i , e ) )
if ( thermalState ( material_homogenizationAt ( e ) ) % sizeState > 0 ) &
2019-03-10 15:53:39 +05:30
thermalState ( material_homogenizationAt ( e ) ) % State ( : , mappingHomogenization ( 1 , i , e ) ) = &
2019-05-16 21:54:54 +05:30
thermalState ( material_homogenizationAt ( e ) ) % subState0 ( : , mappingHomogenization ( 1 , i , e ) )
if ( damageState ( material_homogenizationAt ( e ) ) % sizeState > 0 ) &
2019-03-10 15:53:39 +05:30
damageState ( material_homogenizationAt ( e ) ) % State ( : , mappingHomogenization ( 1 , i , e ) ) = &
2019-05-16 21:54:54 +05:30
damageState ( material_homogenizationAt ( e ) ) % subState0 ( : , mappingHomogenization ( 1 , i , e ) )
2015-10-15 00:06:19 +05:30
endif
2013-01-29 15:58:01 +05:30
endif converged
2015-10-15 00:06:19 +05:30
2012-12-14 20:00:08 +05:30
if ( materialpoint_subStep ( i , e ) > subStepMinHomog ) then
materialpoint_requested ( i , e ) = . true .
2018-09-20 10:57:12 +05:30
materialpoint_subF ( 1 : 3 , 1 : 3 , i , e ) = materialpoint_subF0 ( 1 : 3 , 1 : 3 , i , e ) &
+ materialpoint_subStep ( i , e ) * ( materialpoint_F ( 1 : 3 , 1 : 3 , i , e ) &
- materialpoint_F0 ( 1 : 3 , 1 : 3 , i , e ) )
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
materialpoint_subdt ( i , e ) = materialpoint_subStep ( i , e ) * dt
2013-01-29 15:58:01 +05:30
materialpoint_doneAndHappy ( 1 : 2 , i , e ) = [ . false . , . true . ]
2009-05-07 21:57:36 +05:30
endif
2013-01-29 15:58:01 +05:30
enddo IpLooping1
enddo elementLooping1
2010-11-03 22:52:48 +05:30
!$OMP END PARALLEL DO
2009-05-07 21:57:36 +05:30
2019-04-06 00:15:56 +05:30
NiterationMPstate = 0
2015-10-15 00:06:19 +05:30
2013-01-29 15:58:01 +05:30
convergenceLooping : do while ( . not . terminallyIll . and . &
2010-09-02 02:34:02 +05:30
any ( materialpoint_requested ( : , FEsolving_execELem ( 1 ) : FEsolving_execElem ( 2 ) ) &
2009-05-07 21:57:36 +05:30
. and . . not . materialpoint_doneAndHappy ( 1 , : , FEsolving_execELem ( 1 ) : FEsolving_execElem ( 2 ) ) &
2010-09-02 02:34:02 +05:30
) . and . &
2013-01-29 15:58:01 +05:30
NiterationMPstate < nMPstate )
2009-08-11 22:01:57 +05:30
NiterationMPstate = NiterationMPstate + 1
2009-05-07 21:57:36 +05:30
2013-01-29 15:58:01 +05:30
!--------------------------------------------------------------------------------------------------
! deformation partitioning
2015-10-15 00:06:19 +05:30
! based on materialpoint_subF0,.._subF,crystallite_partionedF0, and homogenization_state,
2009-05-07 21:57:36 +05:30
! results in crystallite_partionedF
2010-11-03 22:52:48 +05:30
!$OMP PARALLEL DO PRIVATE(myNgrains)
2013-01-29 15:58:01 +05:30
elementLooping2 : do e = FEsolving_execElem ( 1 ) , FEsolving_execElem ( 2 )
2019-06-07 00:44:37 +05:30
myNgrains = homogenization_Ngrains ( material_homogenizationAt ( e ) )
2013-01-29 15:58:01 +05:30
IpLooping2 : do i = FEsolving_execIP ( 1 , e ) , FEsolving_execIP ( 2 , e )
if ( materialpoint_requested ( i , e ) . and . & ! process requested but...
. not . materialpoint_doneAndHappy ( 1 , i , e ) ) then ! ...not yet done material points
2019-02-02 15:23:55 +05:30
call partitionDeformation ( i , e ) ! partition deformation onto constituents
2013-01-29 15:58:01 +05:30
crystallite_dt ( 1 : myNgrains , i , e ) = materialpoint_subdt ( i , e ) ! propagate materialpoint dt to grains
crystallite_requested ( 1 : myNgrains , i , e ) = . true . ! request calculation for constituents
2009-08-27 17:40:06 +05:30
else
2013-01-29 15:58:01 +05:30
crystallite_requested ( 1 : myNgrains , i , e ) = . false . ! calculation for constituents not required anymore
2009-05-07 21:57:36 +05:30
endif
2013-01-29 15:58:01 +05:30
enddo IpLooping2
enddo elementLooping2
2010-11-03 22:52:48 +05:30
!$OMP END PARALLEL DO
2015-10-15 00:06:19 +05:30
2013-01-29 15:58:01 +05:30
!--------------------------------------------------------------------------------------------------
! crystallite integration
2009-05-07 21:57:36 +05:30
! based on crystallite_partionedF0,.._partionedF
! incrementing by crystallite_dt
2019-02-02 15:23:55 +05:30
materialpoint_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
2010-11-03 22:52:48 +05:30
!$OMP PARALLEL DO
2013-01-29 15:58:01 +05:30
elementLooping3 : do e = FEsolving_execElem ( 1 ) , FEsolving_execElem ( 2 )
IpLooping3 : do i = FEsolving_execIP ( 1 , e ) , FEsolving_execIP ( 2 , e )
2009-05-07 21:57:36 +05:30
if ( materialpoint_requested ( i , e ) . and . &
. not . materialpoint_doneAndHappy ( 1 , i , e ) ) then
2019-01-18 20:00:50 +05:30
if ( . not . materialpoint_converged ( i , e ) ) then
2013-01-29 15:58:01 +05:30
materialpoint_doneAndHappy ( 1 : 2 , i , e ) = [ . true . , . false . ]
2009-08-11 22:01:57 +05:30
else
2019-01-13 14:18:47 +05:30
materialpoint_doneAndHappy ( 1 : 2 , i , e ) = updateState ( i , e )
2015-05-28 22:32:23 +05:30
materialpoint_converged ( i , e ) = all ( materialpoint_doneAndHappy ( 1 : 2 , i , e ) ) ! converged if done and happy
2009-08-11 22:01:57 +05:30
endif
2009-05-07 21:57:36 +05:30
endif
2013-01-29 15:58:01 +05:30
enddo IpLooping3
enddo elementLooping3
2010-11-03 22:52:48 +05:30
!$OMP END PARALLEL DO
2009-05-07 21:57:36 +05:30
2013-01-29 15:58:01 +05:30
enddo convergenceLooping
2009-05-07 21:57:36 +05:30
2019-04-06 00:15:56 +05:30
NiterationHomog = NiterationHomog + 1
2009-08-11 22:01:57 +05:30
2013-01-29 15:58:01 +05:30
enddo cutBackLooping
2019-01-18 20:00:50 +05:30
if ( updateJaco ) call crystallite_stressTangent
2009-05-07 21:57:36 +05:30
2015-10-15 00:06:19 +05:30
if ( . not . terminallyIll ) then
2013-01-29 15:58:01 +05:30
call crystallite_orientations ( ) ! calculate crystal orientations
2010-11-03 22:52:48 +05:30
!$OMP PARALLEL DO
2013-01-29 15:58:01 +05:30
elementLooping4 : do e = FEsolving_execElem ( 1 ) , FEsolving_execElem ( 2 )
IpLooping4 : do i = FEsolving_execIP ( 1 , e ) , FEsolving_execIP ( 2 , e )
2019-01-13 14:18:47 +05:30
call averageStressAndItsTangent ( i , e )
2013-01-29 15:58:01 +05:30
enddo IpLooping4
enddo elementLooping4
2010-11-03 22:52:48 +05:30
!$OMP END PARALLEL DO
2010-09-02 02:34:02 +05:30
else
2013-01-29 15:58:01 +05:30
write ( 6 , '(/,a,/)' ) '<< HOMOG >> Material Point terminally ill'
2010-03-19 19:44:08 +05:30
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 parallelized calculation of result array at material points
!--------------------------------------------------------------------------------------------------
2013-10-19 00:27:28 +05:30
subroutine materialpoint_postResults
2009-05-07 21:57:36 +05:30
2019-04-06 00:15:56 +05:30
integer :: &
2013-01-29 15:58:01 +05:30
thePos , &
theSize , &
myNgrains , &
myCrystallite , &
g , & !< grain number
i , & !< integration point number
e !< element number
2009-05-07 21:57:36 +05:30
2011-08-01 23:40:55 +05:30
!$OMP PARALLEL DO PRIVATE(myNgrains,myCrystallite,thePos,theSize)
2013-01-29 15:58:01 +05:30
elementLooping : do e = FEsolving_execElem ( 1 ) , FEsolving_execElem ( 2 )
2019-06-07 00:44:37 +05:30
myNgrains = homogenization_Ngrains ( material_homogenizationAt ( e ) )
2019-06-07 02:19:17 +05:30
myCrystallite = microstructure_crystallite ( discretization_microstructureAt ( e ) )
2013-01-29 15:58:01 +05:30
IpLooping : do i = FEsolving_execIP ( 1 , e ) , FEsolving_execIP ( 2 , e )
2019-04-06 00:15:56 +05:30
thePos = 0
2015-10-15 00:06:19 +05:30
2019-03-10 15:53:39 +05:30
theSize = homogState ( material_homogenizationAt ( e ) ) % sizePostResults &
+ thermalState ( material_homogenizationAt ( e ) ) % sizePostResults &
+ damageState ( material_homogenizationAt ( e ) ) % sizePostResults
2012-08-10 21:28:17 +05:30
materialpoint_results ( thePos + 1 , i , e ) = real ( theSize , pReal ) ! tell size of homogenization results
2019-04-06 00:15:56 +05:30
thePos = thePos + 1
2011-03-29 12:57:19 +05:30
2019-04-06 00:15:56 +05:30
if ( theSize > 0 ) then ! any homogenization results to mention?
2019-01-13 14:18:47 +05:30
materialpoint_results ( thePos + 1 : thePos + theSize , i , e ) = postResults ( i , e ) ! tell homogenization results
2011-08-01 23:40:55 +05:30
thePos = thePos + theSize
2009-05-07 21:57:36 +05:30
endif
2014-09-26 16:04:36 +05:30
2012-08-10 21:28:17 +05:30
materialpoint_results ( thePos + 1 , i , e ) = real ( myNgrains , pReal ) ! tell number of grains at materialpoint
2019-04-06 00:15:56 +05:30
thePos = thePos + 1
2011-11-23 14:39:00 +05:30
2013-01-29 15:58:01 +05:30
grainLooping : do g = 1 , myNgrains
2014-09-26 20:46:10 +05:30
theSize = 1 + crystallite_sizePostResults ( myCrystallite ) + &
2019-06-15 17:27:24 +05:30
1 + plasticState ( material_phaseAt ( g , e ) ) % sizePostResults + & !ToDo
sum ( sourceState ( material_phaseAt ( g , e ) ) % p ( : ) % sizePostResults )
2014-03-12 13:03:51 +05:30
materialpoint_results ( thePos + 1 : thePos + theSize , i , e ) = crystallite_postResults ( g , i , e ) ! tell crystallite results
2011-08-01 23:40:55 +05:30
thePos = thePos + theSize
2013-01-29 15:58:01 +05:30
enddo grainLooping
enddo IpLooping
enddo elementLooping
2010-11-03 22:52:48 +05:30
!$OMP END PARALLEL DO
2009-05-07 21:57:36 +05:30
2012-08-10 21:28:17 +05:30
end subroutine materialpoint_postResults
2015-10-15 00:06:19 +05:30
2012-08-10 21:28:17 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief partition material point def grad onto constituents
!--------------------------------------------------------------------------------------------------
2019-01-13 14:18:47 +05:30
subroutine partitionDeformation ( ip , el )
2009-05-07 21:57:36 +05:30
2019-04-06 00:15:56 +05:30
integer , intent ( in ) :: &
2013-10-11 21:31:53 +05:30
ip , & !< integration point
el !< element number
2014-03-14 04:50:50 +05:30
2019-06-07 00:44:37 +05:30
chosenHomogenization : select case ( homogenization_type ( material_homogenizationAt ( el ) ) )
2014-03-14 04:50:50 +05:30
case ( HOMOGENIZATION_NONE_ID ) chosenHomogenization
2019-01-13 13:44:23 +05:30
crystallite_partionedF ( 1 : 3 , 1 : 3 , 1 , ip , el ) = materialpoint_subF ( 1 : 3 , 1 : 3 , ip , el )
2014-03-14 04:50:50 +05:30
2013-11-27 13:34:05 +05:30
case ( HOMOGENIZATION_ISOSTRAIN_ID ) chosenHomogenization
2019-04-06 01:12:07 +05:30
call mech_isostrain_partitionDeformation ( &
2019-06-07 00:44:37 +05:30
crystallite_partionedF ( 1 : 3 , 1 : 3 , 1 : homogenization_Ngrains ( material_homogenizationAt ( el ) ) , ip , el ) , &
2019-01-13 13:44:23 +05:30
materialpoint_subF ( 1 : 3 , 1 : 3 , ip , el ) )
2018-11-03 21:10:17 +05:30
2013-11-27 13:34:05 +05:30
case ( HOMOGENIZATION_RGC_ID ) chosenHomogenization
2019-05-18 10:53:46 +05:30
call mech_RGC_partitionDeformation ( &
2019-06-07 00:44:37 +05:30
crystallite_partionedF ( 1 : 3 , 1 : 3 , 1 : homogenization_Ngrains ( material_homogenizationAt ( el ) ) , ip , el ) , &
2014-08-21 23:18:20 +05:30
materialpoint_subF ( 1 : 3 , 1 : 3 , ip , el ) , &
ip , &
el )
2013-01-29 15:58:01 +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
!--------------------------------------------------------------------------------------------------
2019-01-13 14:18:47 +05:30
function updateState ( ip , el )
2013-11-27 13:34:05 +05:30
2019-04-06 00:15:56 +05:30
integer , intent ( in ) :: &
2013-10-11 21:31:53 +05:30
ip , & !< integration point
el !< element number
2019-01-13 14:18:47 +05:30
logical , dimension ( 2 ) :: updateState
2015-10-15 00:06:19 +05:30
2019-01-13 14:18:47 +05:30
updateState = . true .
2019-06-07 00:44:37 +05:30
chosenHomogenization : select case ( homogenization_type ( material_homogenizationAt ( el ) ) )
2013-11-27 13:34:05 +05:30
case ( HOMOGENIZATION_RGC_ID ) chosenHomogenization
2019-01-13 14:18:47 +05:30
updateState = &
updateState . and . &
2019-06-07 00:44:37 +05:30
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 ) , &
2019-05-18 10:53:46 +05:30
materialpoint_subF ( 1 : 3 , 1 : 3 , ip , el ) , &
materialpoint_subdt ( ip , el ) , &
2019-06-07 00:44:37 +05:30
crystallite_dPdF ( 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 , 1 : homogenization_Ngrains ( material_homogenizationAt ( el ) ) , ip , el ) , &
2019-05-18 10:53:46 +05:30
ip , &
el )
2013-01-29 15:58:01 +05:30
end select chosenHomogenization
2009-05-07 21:57:36 +05:30
2019-06-07 00:44:37 +05:30
chosenThermal : select case ( thermal_type ( material_homogenizationAt ( el ) ) )
2015-05-28 22:32:23 +05:30
case ( THERMAL_adiabatic_ID ) chosenThermal
2019-01-13 14:18:47 +05:30
updateState = &
updateState . and . &
2015-05-28 22:32:23 +05:30
thermal_adiabatic_updateState ( materialpoint_subdt ( ip , el ) , &
ip , &
el )
end select chosenThermal
2019-06-07 00:44:37 +05:30
chosenDamage : select case ( damage_type ( material_homogenizationAt ( el ) ) )
2015-05-28 22:32:23 +05:30
case ( DAMAGE_local_ID ) chosenDamage
2019-01-13 14:18:47 +05:30
updateState = &
updateState . and . &
2015-05-28 22:32:23 +05:30
damage_local_updateState ( materialpoint_subdt ( ip , el ) , &
ip , &
el )
end select chosenDamage
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
2019-04-06 00:15:56 +05:30
integer , intent ( in ) :: &
2013-10-11 21:31:53 +05:30
ip , & !< integration point
el !< element number
2015-10-15 00:06:19 +05:30
2019-06-07 00:44:37 +05:30
chosenHomogenization : select case ( homogenization_type ( material_homogenizationAt ( el ) ) )
2014-03-14 04:50:50 +05:30
case ( HOMOGENIZATION_NONE_ID ) chosenHomogenization
2019-01-13 13:44:23 +05:30
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 )
2015-10-15 00:06:19 +05:30
2013-11-27 13:34:05 +05:30
case ( HOMOGENIZATION_ISOSTRAIN_ID ) chosenHomogenization
2019-04-06 01:12:07 +05:30
call mech_isostrain_averageStressAndItsTangent ( &
2013-10-16 18:34:59 +05:30
materialpoint_P ( 1 : 3 , 1 : 3 , ip , el ) , &
materialpoint_dPdF ( 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 , ip , el ) , &
2019-06-07 00:44:37 +05:30
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 ) ) )
2018-11-03 21:10:17 +05:30
2013-11-27 13:34:05 +05:30
case ( HOMOGENIZATION_RGC_ID ) chosenHomogenization
2019-05-18 10:53:46 +05:30
call mech_RGC_averageStressAndItsTangent ( &
2013-10-16 18:34:59 +05:30
materialpoint_P ( 1 : 3 , 1 : 3 , ip , el ) , &
materialpoint_dPdF ( 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 , ip , el ) , &
2019-06-07 00:44:37 +05:30
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 ) ) )
2013-01-29 15:58:01 +05:30
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
2012-08-10 21:28:17 +05:30
!--------------------------------------------------------------------------------------------------
2015-10-15 00:06:19 +05:30
!> @brief return array of homogenization results for post file inclusion. call only,
2013-01-29 15:58:01 +05:30
!> if homogenization_sizePostResults(i,e) > 0 !!
2012-08-10 21:28:17 +05:30
!--------------------------------------------------------------------------------------------------
2019-01-13 14:18:47 +05:30
function postResults ( ip , el )
2015-10-15 00:06:19 +05:30
2019-04-06 00:15:56 +05:30
integer , intent ( in ) :: &
2013-10-11 21:31:53 +05:30
ip , & !< integration point
el !< element number
2019-03-10 15:53:39 +05:30
real ( pReal ) , dimension ( homogState ( material_homogenizationAt ( el ) ) % sizePostResults &
+ thermalState ( material_homogenizationAt ( el ) ) % sizePostResults &
+ damageState ( material_homogenizationAt ( el ) ) % sizePostResults ) :: &
2019-01-13 14:18:47 +05:30
postResults
2019-04-06 00:15:56 +05:30
integer :: &
2019-01-13 03:37:35 +05:30
startPos , endPos , &
2019-05-17 10:06:30 +05:30
homog
2015-10-15 00:06:19 +05:30
2015-05-28 22:32:23 +05:30
2019-01-13 14:18:47 +05:30
postResults = 0.0_pReal
2019-05-17 10:06:30 +05:30
startPos = 1
endPos = thermalState ( material_homogenizationAt ( el ) ) % sizePostResults
2019-06-07 00:44:37 +05:30
chosenThermal : select case ( thermal_type ( material_homogenizationAt ( el ) ) )
2015-05-28 22:32:23 +05:30
case ( THERMAL_adiabatic_ID ) chosenThermal
2019-03-10 15:53:39 +05:30
homog = material_homogenizationAt ( el )
2019-02-13 03:22:21 +05:30
postResults ( startPos : endPos ) = &
2018-12-31 02:24:50 +05:30
thermal_adiabatic_postResults ( homog , thermal_typeInstance ( homog ) , thermalMapping ( homog ) % p ( ip , el ) )
2015-05-28 22:32:23 +05:30
case ( THERMAL_conduction_ID ) chosenThermal
2019-03-10 15:53:39 +05:30
homog = material_homogenizationAt ( el )
2019-02-13 03:22:21 +05:30
postResults ( startPos : endPos ) = &
2018-12-31 02:24:50 +05:30
thermal_conduction_postResults ( homog , thermal_typeInstance ( homog ) , thermalMapping ( homog ) % p ( ip , el ) )
2019-02-13 03:22:21 +05:30
2015-05-28 22:32:23 +05:30
end select chosenThermal
2019-04-06 00:15:56 +05:30
startPos = endPos + 1
2019-03-10 15:53:39 +05:30
endPos = endPos + damageState ( material_homogenizationAt ( el ) ) % sizePostResults
2019-06-07 00:44:37 +05:30
chosenDamage : select case ( damage_type ( material_homogenizationAt ( el ) ) )
2015-05-28 22:32:23 +05:30
case ( DAMAGE_local_ID ) chosenDamage
2019-01-13 14:18:47 +05:30
postResults ( startPos : endPos ) = damage_local_postResults ( ip , el )
2015-05-28 22:32:23 +05:30
case ( DAMAGE_nonlocal_ID ) chosenDamage
2019-01-13 14:18:47 +05:30
postResults ( startPos : endPos ) = damage_nonlocal_postResults ( ip , el )
2019-01-13 03:37:35 +05:30
2015-05-28 22:32:23 +05:30
end select chosenDamage
2019-01-13 14:18:47 +05:30
end function postResults
2014-09-26 16:04:36 +05:30
2019-04-30 22:15:16 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief writes homogenization results to HDF5 output file
!--------------------------------------------------------------------------------------------------
subroutine homogenization_results
#if defined(PETSc) || defined(DAMASK_HDF5)
use material , only : &
material_homogenization_type = > homogenization_type
integer :: p
character ( len = 256 ) :: group
do p = 1 , size ( config_name_homogenization )
group = trim ( 'current/materialpoint' ) / / '/' / / trim ( config_name_homogenization ( p ) )
call HDF5_closeGroup ( results_addGroup ( group ) )
group = trim ( group ) / / '/mech'
call HDF5_closeGroup ( results_addGroup ( group ) )
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
enddo
#endif
end subroutine homogenization_results
2012-08-10 21:28:17 +05:30
end module homogenization