DAMASK_EICMD/src/material.f90

329 lines
14 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @author Franz Roters, Max-Planck-Institut für Eisenforschung GmbH
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
2018-06-26 22:39:46 +05:30
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
2020-09-20 23:16:33 +05:30
!> @brief Defines phase and homogenization
!--------------------------------------------------------------------------------------------------
module material
2019-06-05 13:32:55 +05:30
use prec
use math
use config
use results
use IO
2019-06-07 18:01:42 +05:30
use rotations
2019-06-07 09:48:42 +05:30
use discretization
2019-09-20 00:10:59 +05:30
implicit none
private
2020-03-16 22:59:15 +05:30
2020-03-17 12:47:14 +05:30
enum, bind(c); enumerator :: &
THERMAL_ISOTHERMAL_ID, &
THERMAL_CONDUCTION_ID, &
DAMAGE_NONE_ID, &
DAMAGE_NONLOCAL_ID, &
HOMOGENIZATION_UNDEFINED_ID, &
HOMOGENIZATION_NONE_ID, &
HOMOGENIZATION_ISOSTRAIN_ID, &
HOMOGENIZATION_RGC_ID
2019-09-20 00:10:59 +05:30
end enum
2020-03-16 22:59:15 +05:30
character(len=:), public, protected, allocatable, dimension(:) :: &
2020-08-15 19:32:10 +05:30
material_name_phase, & !< name of each phase
material_name_homogenization !< name of each homogenization
2019-09-20 00:10:59 +05:30
integer(kind(THERMAL_isothermal_ID)), dimension(:), allocatable, public, protected :: &
thermal_type !< thermal transport model
integer(kind(DAMAGE_none_ID)), dimension(:), allocatable, public, protected :: &
damage_type !< nonlocal damage model
integer(kind(HOMOGENIZATION_undefined_ID)), dimension(:), allocatable, public, protected :: &
homogenization_type !< type of each homogenization
2020-03-16 22:59:15 +05:30
2019-09-20 00:10:59 +05:30
integer, public, protected :: &
homogenization_maxNconstituents !< max number of grains in any USED homogenization
2020-03-16 22:59:15 +05:30
2019-09-20 00:10:59 +05:30
integer, dimension(:), allocatable, public, protected :: &
homogenization_Nconstituents, & !< number of grains in each homogenization
2019-09-20 00:10:59 +05:30
homogenization_typeInstance, & !< instance of particular type of each homogenization
thermal_typeInstance, & !< instance of particular type of each thermal transport
2019-11-25 13:14:44 +05:30
damage_typeInstance !< instance of particular type of each nonlocal damage
2020-03-16 22:59:15 +05:30
2019-09-20 00:10:59 +05:30
real(pReal), dimension(:), allocatable, public, protected :: &
thermal_initialT !< initial temperature per each homogenization
2019-09-20 00:10:59 +05:30
integer, dimension(:), allocatable, public, protected :: & ! (elem)
material_homogenizationAt !< homogenization ID of each element
2020-12-16 00:25:55 +05:30
integer, dimension(:,:), allocatable, public, protected :: & ! (ip,elem)
2019-09-20 00:10:59 +05:30
material_homogenizationMemberAt !< position of the element within its homogenization instance
2020-10-24 18:26:03 +05:30
integer, dimension(:,:), allocatable, public, protected :: & ! (constituent,elem)
2019-09-20 00:10:59 +05:30
material_phaseAt !< phase ID of each element
integer, dimension(:,:,:), allocatable, public, protected :: & ! (constituent,IP,elem)
2019-09-20 00:10:59 +05:30
material_phaseMemberAt !< position of the element within its phase instance
2019-09-20 00:10:59 +05:30
type(tState), allocatable, dimension(:), public :: &
homogState, &
2021-01-08 11:40:38 +05:30
damageState_h
2020-03-16 22:59:15 +05:30
type(Rotation), dimension(:,:,:), allocatable, public, protected :: &
material_orientation0 !< initial orientation of each grain,IP,element
2020-03-16 22:59:15 +05:30
2019-09-20 00:10:59 +05:30
type(group_float), allocatable, dimension(:), public :: &
temperature, & !< temperature field
damage, & !< damage field
temperatureRate !< temperature change rate field
2020-03-16 22:59:15 +05:30
2019-09-20 00:10:59 +05:30
public :: &
material_init, &
2020-03-17 12:47:14 +05:30
THERMAL_ISOTHERMAL_ID, &
THERMAL_CONDUCTION_ID, &
DAMAGE_NONE_ID, &
DAMAGE_NONLOCAL_ID, &
HOMOGENIZATION_NONE_ID, &
HOMOGENIZATION_ISOSTRAIN_ID, &
2019-09-20 00:10:59 +05:30
HOMOGENIZATION_RGC_ID
contains
!--------------------------------------------------------------------------------------------------
!> @brief parses material configuration file
!--------------------------------------------------------------------------------------------------
subroutine material_init(restart)
logical, intent(in) :: restart
2020-09-22 16:39:12 +05:30
print'(/,a)', ' <<<+- material init -+>>>'; flush(IO_STDOUT)
2020-03-16 22:59:15 +05:30
2020-08-15 19:32:10 +05:30
2020-10-01 19:37:50 +05:30
call material_parseMaterial
print*, 'Material parsed'
2020-03-16 22:59:15 +05:30
call material_parseHomogenization
2020-09-19 12:54:27 +05:30
print*, 'Homogenization parsed'
2020-03-16 22:59:15 +05:30
2020-10-24 18:26:03 +05:30
allocate(homogState (size(material_name_homogenization)))
2021-01-08 11:40:38 +05:30
allocate(damageState_h (size(material_name_homogenization)))
2020-03-16 22:59:15 +05:30
2020-10-24 18:26:03 +05:30
allocate(temperature (size(material_name_homogenization)))
allocate(damage (size(material_name_homogenization)))
allocate(temperatureRate (size(material_name_homogenization)))
2019-06-13 03:01:46 +05:30
if (.not. restart) then
call results_openJobFile
call results_mapping_phase(material_phaseAt,material_phaseMemberAt,material_name_phase)
call results_mapping_homogenization(material_homogenizationAt,material_homogenizationMemberAt,material_name_homogenization)
call results_closeJobFile
endif
2019-04-04 16:55:29 +05:30
end subroutine material_init
!--------------------------------------------------------------------------------------------------
2018-06-11 03:46:48 +05:30
!> @brief parses the homogenization part from the material configuration
! ToDo: This should be done in homogenization
!--------------------------------------------------------------------------------------------------
subroutine material_parseHomogenization
2015-10-14 00:22:01 +05:30
2020-08-15 19:32:10 +05:30
class(tNode), pointer :: &
material_homogenization, &
homog, &
homogMech, &
homogThermal, &
homogDamage
2020-01-23 18:30:56 +05:30
2020-08-15 19:32:10 +05:30
integer :: h
2020-03-16 22:59:15 +05:30
material_homogenization => config_material%get('homogenization')
2020-08-15 19:32:10 +05:30
2020-10-24 18:26:03 +05:30
allocate(homogenization_type(size(material_name_homogenization)), source=HOMOGENIZATION_undefined_ID)
allocate(thermal_type(size(material_name_homogenization)), source=THERMAL_isothermal_ID)
allocate(damage_type (size(material_name_homogenization)), source=DAMAGE_none_ID)
allocate(homogenization_typeInstance(size(material_name_homogenization)), source=0)
allocate(thermal_typeInstance(size(material_name_homogenization)), source=0)
allocate(damage_typeInstance(size(material_name_homogenization)), source=0)
allocate(thermal_initialT(size(material_name_homogenization)), source=300.0_pReal)
2020-08-15 19:32:10 +05:30
2020-10-24 18:26:03 +05:30
do h=1, size(material_name_homogenization)
2020-08-15 19:32:10 +05:30
homog => material_homogenization%get(h)
2020-11-18 01:54:40 +05:30
homogMech => homog%get('mechanics')
2020-08-15 19:32:10 +05:30
select case (homogMech%get_asString('type'))
case('none')
2019-09-20 00:10:59 +05:30
homogenization_type(h) = HOMOGENIZATION_NONE_ID
2020-08-15 19:32:10 +05:30
case('isostrain')
2019-09-20 00:10:59 +05:30
homogenization_type(h) = HOMOGENIZATION_ISOSTRAIN_ID
2020-08-15 19:32:10 +05:30
case('RGC')
2019-09-20 00:10:59 +05:30
homogenization_type(h) = HOMOGENIZATION_RGC_ID
case default
2020-08-15 19:32:10 +05:30
call IO_error(500,ext_msg=homogMech%get_asString('type'))
2019-09-20 00:10:59 +05:30
end select
2020-03-16 22:59:15 +05:30
2019-09-20 00:10:59 +05:30
homogenization_typeInstance(h) = count(homogenization_type==homogenization_type(h))
2020-03-16 22:59:15 +05:30
2020-08-15 19:32:10 +05:30
if(homog%contains('thermal')) then
homogThermal => homog%get('thermal')
thermal_initialT(h) = homogThermal%get_asFloat('T_0',defaultVal=300.0_pReal)
select case (homogThermal%get_asString('type'))
case('isothermal')
thermal_type(h) = THERMAL_isothermal_ID
case('conduction')
thermal_type(h) = THERMAL_conduction_ID
case default
call IO_error(500,ext_msg=homogThermal%get_asString('type'))
end select
2019-09-20 00:10:59 +05:30
endif
2020-03-16 22:59:15 +05:30
2020-08-15 19:32:10 +05:30
if(homog%contains('damage')) then
homogDamage => homog%get('damage')
select case (homogDamage%get_asString('type'))
case('none')
damage_type(h) = DAMAGE_none_ID
case('nonlocal')
damage_type(h) = DAMAGE_nonlocal_ID
case default
call IO_error(500,ext_msg=homogDamage%get_asString('type'))
end select
2019-09-20 00:10:59 +05:30
endif
enddo
2020-03-16 22:59:15 +05:30
2020-10-24 18:26:03 +05:30
do h=1, size(material_name_homogenization)
2019-09-20 02:10:03 +05:30
homogenization_typeInstance(h) = count(homogenization_type(1:h) == homogenization_type(h))
thermal_typeInstance(h) = count(thermal_type (1:h) == thermal_type (h))
damage_typeInstance(h) = count(damage_type (1:h) == damage_type (h))
2019-09-20 00:10:59 +05:30
enddo
2020-03-16 22:59:15 +05:30
end subroutine material_parseHomogenization
!--------------------------------------------------------------------------------------------------
2020-10-01 16:13:05 +05:30
!> @brief parses the material part in the material configuration file
!--------------------------------------------------------------------------------------------------
2020-10-01 19:37:50 +05:30
subroutine material_parseMaterial
2015-10-14 00:22:01 +05:30
2020-10-01 16:13:05 +05:30
class(tNode), pointer :: materials, & !> list of materials
material, & !> material definition
constituents, & !> list of constituents
constituent, & !> constituent definition
2020-08-15 19:32:10 +05:30
phases, &
2020-10-24 20:00:48 +05:30
homogenizations, &
homogenization
2020-03-16 22:59:15 +05:30
2020-08-15 19:32:10 +05:30
integer, dimension(:), allocatable :: &
2020-09-19 13:30:49 +05:30
counterPhase, &
counterHomogenization
2020-03-16 22:59:15 +05:30
2020-09-19 13:30:49 +05:30
real(pReal) :: &
frac
2020-08-15 19:32:10 +05:30
integer :: &
2020-10-24 20:00:48 +05:30
e, i, c, &
2020-10-24 21:30:17 +05:30
h
2020-10-24 20:00:48 +05:30
materials => config_material%get('material')
phases => config_material%get('phase')
homogenizations => config_material%get('homogenization')
2020-08-15 19:32:10 +05:30
2020-10-24 21:30:17 +05:30
call sanityCheck(materials, homogenizations)
material_name_phase = getKeys(phases)
material_name_homogenization = getKeys(homogenizations)
2020-10-24 18:18:07 +05:30
allocate(homogenization_Nconstituents(homogenizations%length))
2020-10-24 21:30:17 +05:30
do h=1, homogenizations%length
2020-10-24 20:00:48 +05:30
homogenization => homogenizations%get(h)
homogenization_Nconstituents(h) = homogenization%get_asInt('N_constituents')
2020-10-24 20:00:48 +05:30
enddo
homogenization_maxNconstituents = maxval(homogenization_Nconstituents)
2020-10-24 20:00:48 +05:30
allocate(counterPhase(phases%length),source=0)
allocate(counterHomogenization(homogenizations%length),source=0)
allocate(material_homogenizationAt(discretization_Nelems),source=0)
allocate(material_homogenizationMemberAt(discretization_nIPs,discretization_Nelems),source=0)
allocate(material_phaseAt(homogenization_maxNconstituents,discretization_Nelems),source=0)
allocate(material_phaseMemberAt(homogenization_maxNconstituents,discretization_nIPs,discretization_Nelems),source=0)
allocate(material_orientation0(homogenization_maxNconstituents,discretization_nIPs,discretization_Nelems))
2020-10-24 18:18:07 +05:30
do e = 1, discretization_Nelems
material => materials%get(discretization_materialAt(e))
constituents => material%get('constituents')
material_homogenizationAt(e) = homogenizations%getIndex(material%get_asString('homogenization'))
do i = 1, discretization_nIPs
2020-09-19 13:30:49 +05:30
counterHomogenization(material_homogenizationAt(e)) = counterHomogenization(material_homogenizationAt(e)) + 1
material_homogenizationMemberAt(i,e) = counterHomogenization(material_homogenizationAt(e))
2020-08-15 19:32:10 +05:30
enddo
2020-09-19 13:30:49 +05:30
frac = 0.0_pReal
do c = 1, constituents%length
constituent => constituents%get(c)
frac = frac + constituent%get_asFloat('fraction')
2020-09-19 13:30:49 +05:30
material_phaseAt(c,e) = phases%getIndex(constituent%get_asString('phase'))
do i = 1, discretization_nIPs
2020-09-19 13:30:49 +05:30
counterPhase(material_phaseAt(c,e)) = counterPhase(material_phaseAt(c,e)) + 1
material_phaseMemberAt(c,i,e) = counterPhase(material_phaseAt(c,e))
2020-10-24 18:18:07 +05:30
call material_orientation0(c,i,e)%fromQuaternion(constituent%get_asFloats('O',requiredSize=4)) ! should be done in crystallite
enddo
2020-08-15 19:32:10 +05:30
enddo
2020-09-19 13:30:49 +05:30
if (dNeq(frac,1.0_pReal)) call IO_error(153,ext_msg='constituent')
2020-03-16 22:59:15 +05:30
enddo
2020-10-01 19:37:50 +05:30
end subroutine material_parseMaterial
2018-10-14 23:46:30 +05:30
2020-10-24 21:30:17 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Check if material.yaml is consistent and contains sufficient # of materials
!--------------------------------------------------------------------------------------------------
subroutine sanityCheck(materials,homogenizations)
class(tNode), intent(in) :: materials, &
homogenizations
class(tNode), pointer :: material, &
homogenization, &
constituents
integer :: m
if(maxval(discretization_materialAt) > materials%length) &
call IO_error(155,ext_msg='More materials requested than found in material.yaml')
do m = 1, materials%length
material => materials%get(m)
constituents => material%get('constituents')
homogenization => homogenizations%get(material%get_asString('homogenization'))
if(constituents%length /= homogenization%get_asInt('N_constituents')) call IO_error(148)
enddo
end subroutine sanityCheck
!--------------------------------------------------------------------------------------------------
2020-11-19 02:21:37 +05:30
!> @brief Get all keys from a dictionary
2020-10-24 21:30:17 +05:30
!--------------------------------------------------------------------------------------------------
function getKeys(dict)
class(tNode), intent(in) :: dict
character(len=:), dimension(:), allocatable :: getKeys
character(len=pStringLen), dimension(:), allocatable :: temp
2020-10-24 21:30:17 +05:30
integer :: i,l
2020-10-24 21:30:17 +05:30
allocate(temp(dict%length))
l = 0
2020-10-24 21:30:17 +05:30
do i=1, dict%length
temp(i) = dict%getKey(i)
l = max(len_trim(temp(i)),l)
enddo
2020-12-10 05:01:58 +05:30
allocate(character(l)::getKeys(dict%length))
do i=1, dict%length
getKeys(i) = trim(temp(i))
2020-10-24 21:30:17 +05:30
enddo
end function getKeys
end module material