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
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2012-03-09 01:55:28 +05:30
|
|
|
module homogenization_isostrain
|
2013-01-28 22:06:26 +05:30
|
|
|
use prec, only: &
|
|
|
|
pInt
|
2012-03-09 01:55:28 +05:30
|
|
|
|
2010-02-25 23:09:11 +05:30
|
|
|
implicit none
|
2013-01-28 22:06:26 +05:30
|
|
|
private
|
2013-12-12 22:39:59 +05:30
|
|
|
enum, bind(c)
|
|
|
|
enumerator :: parallel_ID, &
|
|
|
|
average_ID
|
|
|
|
end enum
|
2018-08-25 14:33:43 +05:30
|
|
|
|
2018-08-27 03:41:15 +05:30
|
|
|
type, private :: tParameters !< container type for internal constitutive parameters
|
|
|
|
integer(pInt) :: &
|
|
|
|
Nconstituents
|
|
|
|
integer(kind(average_ID)) :: &
|
|
|
|
mapping
|
|
|
|
end type
|
2013-11-27 17:09:28 +05:30
|
|
|
|
2018-08-27 03:41:15 +05:30
|
|
|
type(tParameters), dimension(:), allocatable, private :: param !< containers of constitutive parameters (len Ninstance)
|
2013-01-28 22:06:26 +05:30
|
|
|
|
|
|
|
public :: &
|
|
|
|
homogenization_isostrain_init, &
|
|
|
|
homogenization_isostrain_partitionDeformation, &
|
2018-08-25 14:33:43 +05:30
|
|
|
homogenization_isostrain_averageStressAndItsTangent
|
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
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-08-25 14:33:43 +05:30
|
|
|
subroutine homogenization_isostrain_init()
|
2018-02-02 17:06:09 +05:30
|
|
|
#if defined(__GFORTRAN__) || __INTEL_COMPILER >= 1800
|
2017-10-05 20:05:34 +05:30
|
|
|
use, intrinsic :: iso_fortran_env, only: &
|
|
|
|
compiler_version, &
|
|
|
|
compiler_options
|
|
|
|
#endif
|
2014-08-21 23:18:20 +05:30
|
|
|
use prec, only: &
|
2014-09-19 23:29:06 +05:30
|
|
|
pReal
|
2014-09-10 19:44:03 +05:30
|
|
|
use debug, only: &
|
|
|
|
debug_HOMOGENIZATION, &
|
|
|
|
debug_level, &
|
|
|
|
debug_levelBasic
|
2018-08-27 03:41:15 +05:30
|
|
|
use IO, only: &
|
|
|
|
IO_timeStamp, &
|
|
|
|
IO_error, &
|
|
|
|
IO_warning
|
|
|
|
use material, only: &
|
|
|
|
homogenization_type, &
|
|
|
|
material_homog, &
|
|
|
|
homogState, &
|
|
|
|
HOMOGENIZATION_ISOSTRAIN_ID, &
|
|
|
|
HOMOGENIZATION_ISOSTRAIN_LABEL, &
|
|
|
|
homogenization_typeInstance
|
|
|
|
use config, only: &
|
|
|
|
config_homogenization
|
2013-10-11 21:31:53 +05:30
|
|
|
|
|
|
|
implicit none
|
|
|
|
integer(pInt) :: &
|
2018-08-25 14:33:43 +05:30
|
|
|
h
|
2013-10-11 21:31:53 +05:30
|
|
|
integer :: &
|
2018-11-03 21:20:43 +05:30
|
|
|
Ninstance
|
2014-09-19 23:29:06 +05:30
|
|
|
integer :: &
|
|
|
|
NofMyHomog ! no pInt (stores a system dependen value from 'count'
|
2013-07-01 12:10:09 +05:30
|
|
|
character(len=65536) :: &
|
2018-08-25 14:33:43 +05:30
|
|
|
tag = ''
|
2018-08-27 03:41:15 +05:30
|
|
|
type(tParameters) :: prm
|
2010-02-25 23:09:11 +05:30
|
|
|
|
2017-11-21 19:38:45 +05:30
|
|
|
write(6,'(/,a)') ' <<<+- homogenization_'//HOMOGENIZATION_ISOSTRAIN_label//' init -+>>>'
|
|
|
|
write(6,'(a15,a)') ' Current time: ',IO_timeStamp()
|
2012-02-01 00:48:55 +05:30
|
|
|
#include "compilation_info.f90"
|
2013-01-09 03:41:59 +05:30
|
|
|
|
2018-11-03 21:20:43 +05:30
|
|
|
Ninstance = count(homogenization_type == HOMOGENIZATION_ISOSTRAIN_ID)
|
|
|
|
if (Ninstance == 0) return
|
2014-09-19 23:29:06 +05:30
|
|
|
|
|
|
|
if (iand(debug_level(debug_HOMOGENIZATION),debug_levelBasic) /= 0_pInt) &
|
2018-11-03 21:20:43 +05:30
|
|
|
write(6,'(a16,1x,i5,/)') '# instances:',Ninstance
|
2014-09-26 16:04:36 +05:30
|
|
|
|
2018-11-03 21:20:43 +05:30
|
|
|
allocate(param(Ninstance)) ! one container of parameters per instance
|
2018-08-25 14:33:43 +05:30
|
|
|
|
|
|
|
do h = 1_pInt, size(homogenization_type)
|
|
|
|
if (homogenization_type(h) /= HOMOGENIZATION_ISOSTRAIN_ID) cycle
|
2018-11-03 21:20:43 +05:30
|
|
|
associate(prm => param(homogenization_typeInstance(h)))
|
2018-08-25 14:33:43 +05:30
|
|
|
|
2018-08-27 03:41:15 +05:30
|
|
|
prm%Nconstituents = config_homogenization(h)%getInt('nconstituents')
|
2018-08-25 14:33:43 +05:30
|
|
|
tag = 'sum'
|
|
|
|
tag = config_homogenization(h)%getString('mapping',defaultVal = tag)
|
|
|
|
select case(trim(tag))
|
2018-08-27 03:41:15 +05:30
|
|
|
case ('sum')
|
|
|
|
prm%mapping = parallel_ID
|
|
|
|
case ('avg')
|
|
|
|
prm%mapping = average_ID
|
2018-08-25 14:33:43 +05:30
|
|
|
case default
|
|
|
|
call IO_error(211_pInt,ext_msg=trim(tag)//' ('//HOMOGENIZATION_isostrain_label//')')
|
|
|
|
end select
|
|
|
|
|
|
|
|
NofMyHomog = count(material_homog == h)
|
|
|
|
|
|
|
|
homogState(h)%sizeState = 0_pInt
|
|
|
|
homogState(h)%sizePostResults = 0_pInt
|
|
|
|
allocate(homogState(h)%state0 (0_pInt,NofMyHomog), source=0.0_pReal)
|
|
|
|
allocate(homogState(h)%subState0(0_pInt,NofMyHomog), source=0.0_pReal)
|
|
|
|
allocate(homogState(h)%state (0_pInt,NofMyHomog), source=0.0_pReal)
|
2018-08-27 03:41:15 +05:30
|
|
|
end associate
|
2014-09-26 16:04:36 +05:30
|
|
|
|
2018-08-25 14:33:43 +05:30
|
|
|
enddo
|
2014-09-10 19:44:03 +05:30
|
|
|
|
2012-03-09 01:55:28 +05:30
|
|
|
end subroutine homogenization_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
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-11-03 21:10:17 +05:30
|
|
|
subroutine homogenization_isostrain_partitionDeformation(F,avgF,instance)
|
2013-01-28 22:06:26 +05:30
|
|
|
use prec, only: &
|
2013-10-16 18:34:59 +05:30
|
|
|
pReal
|
2013-01-28 22:06:26 +05:30
|
|
|
use material, only: &
|
2018-11-03 21:10:17 +05:30
|
|
|
homogenization_maxNgrains
|
2013-01-28 22:06:26 +05:30
|
|
|
|
2010-02-25 23:09:11 +05:30
|
|
|
implicit none
|
2013-10-16 18:34:59 +05:30
|
|
|
real(pReal), dimension (3,3,homogenization_maxNgrains), intent(out) :: F !< partioned def grad per grain
|
|
|
|
real(pReal), dimension (3,3), intent(in) :: avgF !< my average def grad
|
2018-11-03 21:10:17 +05:30
|
|
|
integer(pInt), intent(in) :: instance
|
2018-08-27 03:41:15 +05:30
|
|
|
type(tParameters) :: &
|
|
|
|
prm
|
|
|
|
|
|
|
|
associate(prm => param(instance))
|
|
|
|
F(1:3,1:3,1:prm%Nconstituents) = spread(avgF,3,prm%Nconstituents)
|
2018-08-28 09:57:38 +05:30
|
|
|
if (homogenization_maxNgrains > prm%Nconstituents) &
|
|
|
|
F(1:3,1:3,prm%Nconstituents+1_pInt:homogenization_maxNgrains) = 0.0_pReal
|
2018-08-27 03:41:15 +05:30
|
|
|
end associate
|
2013-10-11 21:31:53 +05:30
|
|
|
|
|
|
|
end subroutine homogenization_isostrain_partitionDeformation
|
2010-02-25 23:09:11 +05:30
|
|
|
|
|
|
|
|
2013-01-28 22:06:26 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief derive average stress and stiffness from constituent quantities
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-11-03 21:10:17 +05:30
|
|
|
subroutine homogenization_isostrain_averageStressAndItsTangent(avgP,dAvgPdAvgF,P,dPdF,instance)
|
2013-01-28 22:06:26 +05:30
|
|
|
use prec, only: &
|
|
|
|
pReal
|
2013-10-11 21:31:53 +05:30
|
|
|
use material, only: &
|
2018-11-03 21:10:17 +05:30
|
|
|
homogenization_maxNgrains
|
2012-03-09 01:55:28 +05:30
|
|
|
|
2010-02-25 23:09:11 +05:30
|
|
|
implicit none
|
2013-10-11 21:31:53 +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
|
|
|
|
real(pReal), dimension (3,3,homogenization_maxNgrains), intent(in) :: P !< array of current grain stresses
|
|
|
|
real(pReal), dimension (3,3,3,3,homogenization_maxNgrains), intent(in) :: dPdF !< array of current grain stiffnesses
|
2018-11-03 21:10:17 +05:30
|
|
|
integer(pInt), intent(in) :: instance
|
2018-08-27 03:41:15 +05:30
|
|
|
type(tParameters) :: &
|
|
|
|
prm
|
2013-04-30 15:10:06 +05:30
|
|
|
|
2018-08-27 03:41:15 +05:30
|
|
|
associate(prm => param(instance))
|
|
|
|
select case (prm%mapping)
|
2013-12-12 22:39:59 +05:30
|
|
|
case (parallel_ID)
|
2013-04-30 15:10:06 +05:30
|
|
|
avgP = sum(P,3)
|
|
|
|
dAvgPdAvgF = sum(dPdF,5)
|
2013-12-12 22:39:59 +05:30
|
|
|
case (average_ID)
|
2018-08-27 03:41:15 +05:30
|
|
|
avgP = sum(P,3) /real(prm%Nconstituents,pReal)
|
|
|
|
dAvgPdAvgF = sum(dPdF,5)/real(prm%Nconstituents,pReal)
|
2013-04-30 15:10:06 +05:30
|
|
|
end select
|
2018-08-27 03:41:15 +05:30
|
|
|
end associate
|
2010-02-25 23:09:11 +05:30
|
|
|
|
2012-03-09 01:55:28 +05:30
|
|
|
end subroutine homogenization_isostrain_averageStressAndItsTangent
|
2010-02-25 23:09:11 +05:30
|
|
|
|
2012-03-09 01:55:28 +05:30
|
|
|
end module homogenization_isostrain
|