DAMASK_EICMD/src/material.f90

236 lines
9.1 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 config
2023-01-19 22:07:45 +05:30
use result
2021-07-24 18:33:26 +05:30
use math
2019-06-05 13:32:55 +05:30
use IO
2019-06-07 18:01:42 +05:30
use rotations
2019-06-07 09:48:42 +05:30
use discretization
2021-02-12 20:01:43 +05:30
use YAML_types
implicit none(type,external)
2019-09-20 00:10:59 +05:30
private
2020-03-16 22:59:15 +05:30
type, public :: tRotationContainer
type(tRotation), dimension(:), allocatable :: data
end type tRotationContainer
type, public :: tTensorContainer
real(pREAL), dimension(:,:,:), allocatable :: data
end type tTensorContainer
2021-07-24 18:33:26 +05:30
2021-05-22 20:51:07 +05:30
type(tRotationContainer), dimension(:), allocatable, public, protected :: material_O_0
type(tTensorContainer), dimension(:), allocatable, public, protected :: material_V_e_0
2021-05-22 20:51:07 +05:30
2021-02-12 20:01:43 +05:30
integer, dimension(:), allocatable, public, protected :: &
2021-03-01 10:46:16 +05:30
homogenization_Nconstituents !< number of grains in each homogenization
2021-05-22 13:03:58 +05:30
integer, public, protected :: &
homogenization_maxNconstituents !< max number of grains in any homogenization
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
integer, dimension(:), allocatable, public, protected :: & ! (cell)
2023-01-23 13:01:59 +05:30
material_ID_homogenization, & !< Number of the homogenization
material_entry_homogenization !< Position in array of used homogenization
integer, dimension(:,:), allocatable, public, protected :: & ! (constituent,cell)
2023-01-23 13:01:59 +05:30
material_ID_phase, & !< Number of the phase
material_entry_phase !< Position in array of used phase
real(pREAL), dimension(:,:), allocatable, public, protected :: &
material_v ! fraction
2019-09-20 00:10:59 +05:30
public :: &
2023-02-28 12:25:34 +05:30
material_init
contains
!--------------------------------------------------------------------------------------------------
!> @brief Parse material configuration file (material.yaml).
!--------------------------------------------------------------------------------------------------
subroutine material_init(restart)
logical, intent(in) :: restart
2021-02-12 20:01:43 +05:30
print'(/,1x,a)', '<<<+- material init -+>>>'; flush(IO_STDOUT)
2020-03-16 22:59:15 +05:30
2020-08-15 19:32:10 +05:30
2022-01-31 19:35:15 +05:30
call parse()
print'(/,1x,a)', 'parsed material.yaml'
2020-03-16 22:59:15 +05:30
if (.not. restart) then
call result_openJobFile()
2023-01-23 13:01:59 +05:30
call result_mapping_phase(material_ID_phase,material_entry_phase,material_name_phase)
call result_mapping_homogenization(material_ID_homogenization,material_entry_homogenization,material_name_homogenization)
call result_closeJobFile()
end if
2019-04-04 16:55:29 +05:30
end subroutine material_init
!--------------------------------------------------------------------------------------------------
!> @brief Parse material.yaml to get the global structure.
!--------------------------------------------------------------------------------------------------
2021-05-22 13:03:58 +05:30
subroutine parse()
2015-10-14 00:22:01 +05:30
type(tList), pointer :: materials, & !> all materials
constituents !> all constituents of a material
type(tDict), pointer :: phases, & !> all phases
homogenizations, & !> all homogenizations
material, & !> material definition
constituent, & !> constituent definition
homogenization
2020-03-16 22:59:15 +05:30
class(tItem), pointer :: item
2020-08-15 19:32:10 +05:30
integer, dimension(:), allocatable :: &
2020-09-19 13:30:49 +05:30
counterPhase, &
counterHomogenization, &
2022-04-19 22:41:40 +05:30
ho_of
integer, dimension(:,:), allocatable :: ph_of
real(pREAL), dimension(:,:), allocatable :: v_of
2020-03-16 22:59:15 +05:30
real(pREAL) :: v
2020-08-15 19:32:10 +05:30
integer :: &
el, ip, &
ho, ph, &
co, ce, &
ma
2020-10-24 20:00:48 +05:30
materials => config_material%get_list('material')
phases => config_material%get_dict('phase')
homogenizations => config_material%get_dict('homogenization')
2020-08-15 19:32:10 +05:30
2022-04-14 23:44:08 +05:30
if (maxval(discretization_materialAt) > materials%length) &
call IO_error(155,ext_msg='More materials requested than found in material.yaml')
material_name_phase = phases%keys()
material_name_homogenization = homogenizations%keys()
2020-10-24 18:18:07 +05:30
allocate(homogenization_Nconstituents(homogenizations%length))
do ho=1, homogenizations%length
homogenization => homogenizations%get_dict(ho)
homogenization_Nconstituents(ho) = homogenization%get_asInt('N_constituents')
end do
homogenization_maxNconstituents = maxval(homogenization_Nconstituents)
2020-10-24 20:00:48 +05:30
allocate(material_v(homogenization_maxNconstituents,discretization_Ncells),source=0.0_pREAL)
2021-01-17 19:08:12 +05:30
2021-07-24 18:33:26 +05:30
allocate(material_O_0(materials%length))
allocate(material_V_e_0(materials%length))
2021-05-22 20:51:07 +05:30
2022-04-19 22:41:40 +05:30
allocate(ho_of(materials%length))
allocate(ph_of(materials%length,homogenization_maxNconstituents),source=-1)
allocate( v_of(materials%length,homogenization_maxNconstituents),source=0.0_pREAL)
2022-04-14 23:44:08 +05:30
2022-12-25 21:06:01 +05:30
! Parse YAML structure. Manual loop over linked list to have O(n) instead of O(n^2) complexity
item => materials%first
do ma = 1, materials%length
material => item%node%asDict()
2023-06-04 10:47:38 +05:30
ho_of(ma) = homogenizations%index(material%get_asStr('homogenization'))
constituents => material%get_list('constituents')
homogenization => homogenizations%get_dict(ho_of(ma))
if (constituents%length /= homogenization%get_asInt('N_constituents')) call IO_error(148)
allocate(material_O_0(ma)%data(constituents%length))
allocate(material_V_e_0(ma)%data(1:3,1:3,constituents%length))
do co = 1, constituents%length
constituent => constituents%get_dict(co)
v_of(ma,co) = constituent%get_asReal('v')
2023-06-04 10:47:38 +05:30
ph_of(ma,co) = phases%index(constituent%get_asStr('phase'))
call material_O_0(ma)%data(co)%fromQuaternion(constituent%get_as1dReal('O',requiredSize=4))
material_V_e_0(ma)%data(1:3,1:3,co) = constituent%get_as2dReal('V_e',defaultVal=math_I3,requiredShape=[3,3])
if (any(dNeq(material_V_e_0(ma)%data(1:3,1:3,co),transpose(material_V_e_0(ma)%data(1:3,1:3,co))))) &
call IO_error(147)
end do
if (dNeq(sum(v_of(ma,:)),1.0_pREAL,1.e-9_pREAL)) call IO_error(153,ext_msg='constituent')
2018-10-14 23:46:30 +05:30
item => item%next
end do
2022-04-14 23:44:08 +05:30
allocate(counterPhase(phases%length),source=0)
allocate(counterHomogenization(homogenizations%length),source=0)
2023-01-23 13:01:59 +05:30
allocate(material_ID_homogenization(discretization_Ncells),source=0)
allocate(material_entry_homogenization(discretization_Ncells),source=0)
2022-04-14 23:44:08 +05:30
2023-01-23 13:01:59 +05:30
allocate(material_ID_phase(homogenization_maxNconstituents,discretization_Ncells),source=0)
allocate(material_entry_phase(homogenization_maxNconstituents,discretization_Ncells),source=0)
2022-04-14 23:44:08 +05:30
! build mappings
do el = 1, discretization_Nelems
2020-10-24 21:30:17 +05:30
ma = discretization_materialAt(el)
2022-04-19 22:41:40 +05:30
ho = ho_of(ma)
2020-10-24 21:30:17 +05:30
do ip = 1, discretization_nIPs
ce = (el-1)*discretization_nIPs + ip
2023-01-23 13:01:59 +05:30
material_ID_homogenization(ce) = ho
counterHomogenization(ho) = counterHomogenization(ho) + 1
2023-01-23 13:01:59 +05:30
material_entry_homogenization(ce) = counterHomogenization(ho)
end do
2020-10-24 21:30:17 +05:30
2022-04-19 22:41:40 +05:30
do co = 1, size(ph_of(ma,:)>0)
2022-04-19 22:41:40 +05:30
v = v_of(ma,co)
ph = ph_of(ma,co)
2020-10-24 21:30:17 +05:30
do ip = 1, discretization_nIPs
ce = (el-1)*discretization_nIPs + ip
2023-01-23 13:01:59 +05:30
material_ID_phase(co,ce) = ph
counterPhase(ph) = counterPhase(ph) + 1
2023-01-23 13:01:59 +05:30
material_entry_phase(co,ce) = counterPhase(ph)
material_v(co,ce) = v
end do
2020-10-24 21:30:17 +05:30
end do
end do
end subroutine parse
2020-10-24 21:30:17 +05:30
2022-04-14 23:44:08 +05:30
#if defined (__GFORTRAN__)
2020-10-24 21:30:17 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief %keys() is broken on gfortran
2020-10-24 21:30:17 +05:30
!--------------------------------------------------------------------------------------------------
function getKeys(dict)
type(tDict), intent(in) :: dict
2023-06-04 10:47:38 +05:30
character(len=:), dimension(:), allocatable :: getKeys
character(len=pSTRLEN), 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%key(i)
l = max(len_trim(temp(i)),l)
end do
2020-12-10 05:01:58 +05:30
allocate(character(l)::getKeys(dict%length))
do i=1, dict%length
getKeys(i) = trim(temp(i))
end do
2020-10-24 21:30:17 +05:30
end function getKeys
#endif
2020-10-24 21:30:17 +05:30
end module material