DAMASK_EICMD/src/material.f90

240 lines
9.3 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
use results
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
2019-09-20 00:10:59 +05:30
implicit none
private
2020-03-16 22:59:15 +05:30
type, public :: tRotationContainer
type(tRotation), dimension(:), allocatable :: data
end type tRotationContainer
type, public :: tTensorContainer
2021-07-24 18:33:26 +05:30
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_F_i_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)
material_homogenizationID, & ! TODO: rename to material_ID_homogenization
material_homogenizationEntry ! TODO: rename to material_entry_homogenization
integer, dimension(:,:), allocatable, public, protected :: & ! (constituent,cell)
material_phaseID, & ! TODO: rename to material_ID_phase
material_phaseEntry ! TODO: rename to material_entry_phase
real(pReal), dimension(:,:), allocatable, public, protected :: &
material_v ! fraction
2019-09-20 00:10:59 +05:30
public :: &
2021-02-12 20:01:43 +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 results_openJobFile
2021-05-22 17:06:23 +05:30
call results_mapping_phase(material_phaseID,material_phaseEntry,material_name_phase)
call results_mapping_homogenization(material_homogenizationID,material_homogenizationEntry,material_name_homogenization)
call results_closeJobFile
end if
2019-04-04 16:55:29 +05:30
end subroutine material_init
!--------------------------------------------------------------------------------------------------
2021-05-22 13:03:58 +05:30
!> @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
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
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('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)
#if defined (__GFORTRAN__)
2020-10-24 21:30:17 +05:30
material_name_phase = getKeys(phases)
material_name_homogenization = getKeys(homogenizations)
#else
material_name_phase = phases%Keys()
material_name_homogenization = homogenizations%Keys()
#endif
2020-10-24 18:18:07 +05:30
allocate(homogenization_Nconstituents(homogenizations%length))
do ho=1, homogenizations%length
homogenization => homogenizations%get(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(counterPhase(phases%length),source=0)
allocate(counterHomogenization(homogenizations%length),source=0)
allocate(material_homogenizationID(discretization_Ncells),source=0)
allocate(material_homogenizationEntry(discretization_Ncells),source=0)
allocate(material_phaseID(homogenization_maxNconstituents,discretization_Ncells),source=0)
allocate(material_phaseEntry(homogenization_maxNconstituents,discretization_Ncells),source=0)
allocate(material_v(homogenization_maxNconstituents,discretization_Ncells),source=0.0_pReal)
2021-01-17 19:08:12 +05:30
do el = 1, discretization_Nelems
material => materials%get(discretization_materialAt(el))
ho = homogenizations%getIndex(material%get_asString('homogenization'))
2021-01-17 19:08:12 +05:30
do ip = 1, discretization_nIPs
ce = (el-1)*discretization_nIPs + ip
material_homogenizationID(ce) = ho
counterHomogenization(ho) = counterHomogenization(ho) + 1
material_homogenizationEntry(ce) = counterHomogenization(ho)
end do
constituents => material%get('constituents')
2021-01-17 19:08:12 +05:30
do co = 1, constituents%length
constituent => constituents%get(co)
v = constituent%get_asFloat('v')
ph = phases%getIndex(constituent%get_asString('phase'))
2021-01-17 19:08:12 +05:30
do ip = 1, discretization_nIPs
ce = (el-1)*discretization_nIPs + ip
material_phaseID(co,ce) = ph
counterPhase(ph) = counterPhase(ph) + 1
material_phaseEntry(co,ce) = counterPhase(ph)
material_v(co,ce) = v
end do
end do
if (dNeq(sum(material_v(1:constituents%length,ce)),1.0_pReal,1.e-9_pReal)) &
call IO_error(153,ext_msg='constituent')
end do
2020-03-16 22:59:15 +05:30
2021-07-24 18:33:26 +05:30
allocate(material_O_0(materials%length))
allocate(material_F_i_0(materials%length))
2021-05-22 20:51:07 +05:30
do ma = 1, materials%length
material => materials%get(ma)
constituents => material%get('constituents')
2021-07-24 18:33:26 +05:30
allocate(material_O_0(ma)%data(constituents%length))
allocate(material_F_i_0(ma)%data(1:3,1:3,constituents%length))
2021-05-22 20:51:07 +05:30
do co = 1, constituents%length
2021-05-22 22:12:18 +05:30
constituent => constituents%get(co)
2021-07-24 18:33:26 +05:30
call material_O_0(ma)%data(co)%fromQuaternion(constituent%get_as1dFloat('O',requiredSize=4))
2021-08-06 00:02:58 +05:30
material_F_i_0(ma)%data(1:3,1:3,co) = constituent%get_as2dFloat('F_i',defaultVal=math_I3,requiredShape=[3,3])
2021-05-22 20:51:07 +05:30
enddo
enddo
2021-05-22 13:03:58 +05:30
end subroutine parse
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) &
2020-10-24 21:30:17 +05:30
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)
end do
2020-10-24 21:30:17 +05:30
end subroutine sanityCheck
#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)
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)
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