2013-01-28 22:06:26 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-11-03 21:10:17 +05:30
|
|
|
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
|
2013-01-28 22:06:26 +05:30
|
|
|
!> @author Franz Roters, Max-Planck-Institut für Eisenforschung GmbH
|
|
|
|
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
|
|
|
|
!> @brief Isostrain (full constraint Taylor assuption) homogenization scheme
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-04-06 01:12:07 +05:30
|
|
|
submodule(homogenization) homogenization_mech_isostrain
|
|
|
|
|
2020-03-17 12:47:14 +05:30
|
|
|
enum, bind(c); enumerator :: &
|
|
|
|
parallel_ID, &
|
|
|
|
average_ID
|
2019-04-06 00:18:20 +05:30
|
|
|
end enum
|
2020-09-13 14:09:17 +05:30
|
|
|
|
2019-04-06 01:12:07 +05:30
|
|
|
type :: tParameters !< container type for internal constitutive parameters
|
2019-04-06 00:18:20 +05:30
|
|
|
integer :: &
|
2020-09-23 05:03:19 +05:30
|
|
|
N_constituents
|
2019-04-06 00:18:20 +05:30
|
|
|
integer(kind(average_ID)) :: &
|
|
|
|
mapping
|
|
|
|
end type
|
2020-09-13 14:09:17 +05:30
|
|
|
|
2019-04-06 01:12:07 +05:30
|
|
|
type(tParameters), dimension(:), allocatable :: param !< containers of constitutive parameters (len Ninstance)
|
2020-09-13 14:09:17 +05:30
|
|
|
|
2010-02-25 23:09:11 +05:30
|
|
|
|
2012-03-09 01:55:28 +05:30
|
|
|
contains
|
2013-01-28 22:06:26 +05:30
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief allocates all neccessary fields, reads information from material configuration file
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-07-02 01:50:22 +05:30
|
|
|
module subroutine mech_isostrain_init
|
2019-05-15 02:42:32 +05:30
|
|
|
|
2019-04-06 00:18:20 +05:30
|
|
|
integer :: &
|
|
|
|
Ninstance, &
|
|
|
|
h, &
|
|
|
|
NofMyHomog
|
2020-08-15 19:32:10 +05:30
|
|
|
class(tNode), pointer :: &
|
|
|
|
material_homogenization, &
|
|
|
|
homog, &
|
|
|
|
homogMech
|
2020-09-13 14:09:17 +05:30
|
|
|
|
2020-09-18 02:27:56 +05:30
|
|
|
print'(/,a)', ' <<<+- homogenization_mech_isostrain init -+>>>'
|
2020-09-13 14:09:17 +05:30
|
|
|
|
2019-04-06 00:18:20 +05:30
|
|
|
Ninstance = count(homogenization_type == HOMOGENIZATION_ISOSTRAIN_ID)
|
2020-09-22 16:39:12 +05:30
|
|
|
print'(a,i2)', ' # instances: ',Ninstance; flush(IO_STDOUT)
|
2020-09-13 14:09:17 +05:30
|
|
|
|
2019-04-06 00:18:20 +05:30
|
|
|
allocate(param(Ninstance)) ! one container of parameters per instance
|
2020-08-15 19:32:10 +05:30
|
|
|
|
2020-09-13 14:09:17 +05:30
|
|
|
material_homogenization => config_material%get('homogenization')
|
2019-04-06 00:18:20 +05:30
|
|
|
do h = 1, size(homogenization_type)
|
|
|
|
if (homogenization_type(h) /= HOMOGENIZATION_ISOSTRAIN_ID) cycle
|
2020-08-15 19:32:10 +05:30
|
|
|
homog => material_homogenization%get(h)
|
|
|
|
homogMech => homog%get('mech')
|
|
|
|
associate(prm => param(homogenization_typeInstance(h)))
|
2020-09-13 14:09:17 +05:30
|
|
|
|
2020-10-07 21:14:54 +05:30
|
|
|
prm%N_constituents = homogenization_Ngrains(h)
|
2020-08-15 19:32:10 +05:30
|
|
|
select case(homogMech%get_asString('mapping',defaultVal = 'sum'))
|
2019-04-06 00:18:20 +05:30
|
|
|
case ('sum')
|
|
|
|
prm%mapping = parallel_ID
|
|
|
|
case ('avg')
|
|
|
|
prm%mapping = average_ID
|
|
|
|
case default
|
2020-08-15 19:32:10 +05:30
|
|
|
call IO_error(211,ext_msg='sum'//' (mech_isostrain)')
|
2019-04-06 00:18:20 +05:30
|
|
|
end select
|
2020-09-13 14:09:17 +05:30
|
|
|
|
2019-04-06 00:18:20 +05:30
|
|
|
NofMyHomog = count(material_homogenizationAt == h)
|
|
|
|
homogState(h)%sizeState = 0
|
|
|
|
allocate(homogState(h)%state0 (0,NofMyHomog))
|
|
|
|
allocate(homogState(h)%subState0(0,NofMyHomog))
|
|
|
|
allocate(homogState(h)%state (0,NofMyHomog))
|
2020-09-13 14:09:17 +05:30
|
|
|
|
2019-04-06 00:18:20 +05:30
|
|
|
end associate
|
2020-09-13 14:09:17 +05:30
|
|
|
|
2019-04-06 00:18:20 +05:30
|
|
|
enddo
|
2014-09-10 19:44:03 +05:30
|
|
|
|
2019-04-06 01:12:07 +05:30
|
|
|
end subroutine mech_isostrain_init
|
2010-02-25 23:09:11 +05:30
|
|
|
|
|
|
|
|
2013-01-28 22:06:26 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief partitions the deformation gradient onto the constituents
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-04-06 01:12:07 +05:30
|
|
|
module subroutine mech_isostrain_partitionDeformation(F,avgF)
|
2020-09-13 14:09:17 +05:30
|
|
|
|
2019-04-06 00:18:20 +05:30
|
|
|
real(pReal), dimension (:,:,:), intent(out) :: F !< partitioned deformation gradient
|
2020-09-13 14:09:17 +05:30
|
|
|
|
2019-04-06 00:18:20 +05:30
|
|
|
real(pReal), dimension (3,3), intent(in) :: avgF !< average deformation gradient at material point
|
2020-09-13 14:09:17 +05:30
|
|
|
|
2019-04-06 00:18:20 +05:30
|
|
|
F = spread(avgF,3,size(F,3))
|
2013-10-11 21:31:53 +05:30
|
|
|
|
2019-04-06 01:12:07 +05:30
|
|
|
end subroutine mech_isostrain_partitionDeformation
|
2010-02-25 23:09:11 +05:30
|
|
|
|
|
|
|
|
2013-01-28 22:06:26 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2020-09-13 14:09:17 +05:30
|
|
|
!> @brief derive average stress and stiffness from constituent quantities
|
2013-01-28 22:06:26 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-04-06 01:12:07 +05:30
|
|
|
module subroutine mech_isostrain_averageStressAndItsTangent(avgP,dAvgPdAvgF,P,dPdF,instance)
|
2020-09-13 14:09:17 +05:30
|
|
|
|
2019-04-06 00:18:20 +05:30
|
|
|
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-09-13 14:09:17 +05:30
|
|
|
|
2019-04-06 00:18:20 +05:30
|
|
|
real(pReal), dimension (:,:,:), intent(in) :: P !< partitioned stresses
|
|
|
|
real(pReal), dimension (:,:,:,:,:), intent(in) :: dPdF !< partitioned stiffnesses
|
2020-09-13 14:09:17 +05:30
|
|
|
integer, intent(in) :: instance
|
|
|
|
|
2019-04-06 00:18:20 +05:30
|
|
|
associate(prm => param(instance))
|
2020-09-13 14:09:17 +05:30
|
|
|
|
2019-04-06 00:18:20 +05:30
|
|
|
select case (prm%mapping)
|
|
|
|
case (parallel_ID)
|
|
|
|
avgP = sum(P,3)
|
|
|
|
dAvgPdAvgF = sum(dPdF,5)
|
|
|
|
case (average_ID)
|
2020-09-23 05:03:19 +05:30
|
|
|
avgP = sum(P,3) /real(prm%N_constituents,pReal)
|
|
|
|
dAvgPdAvgF = sum(dPdF,5)/real(prm%N_constituents,pReal)
|
2019-04-06 00:18:20 +05:30
|
|
|
end select
|
2020-09-13 14:09:17 +05:30
|
|
|
|
2019-04-06 00:18:20 +05:30
|
|
|
end associate
|
2010-02-25 23:09:11 +05:30
|
|
|
|
2019-04-06 01:12:07 +05:30
|
|
|
end subroutine mech_isostrain_averageStressAndItsTangent
|
2010-02-25 23:09:11 +05:30
|
|
|
|
2019-04-06 01:12:07 +05:30
|
|
|
end submodule homogenization_mech_isostrain
|