DAMASK_EICMD/src/geometry_plastic_nonlocal.f90

149 lines
6.4 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
!> @author Christoph Koords, Max-Planck-Institut für Eisenforschung GmbH
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
!> @brief Geometric information about the IP cells needed for the nonlocal
! plasticity model
!--------------------------------------------------------------------------------------------------
module geometry_plastic_nonlocal
use prec
use results
implicit none
2020-01-03 00:17:48 +05:30
public
2020-01-03 00:17:48 +05:30
integer, protected :: &
geometry_plastic_nonlocal_nIPneighbors
2020-01-03 00:17:48 +05:30
integer, dimension(:,:,:,:), allocatable, protected :: &
geometry_plastic_nonlocal_IPneighborhood !< 6 or less neighboring IPs as [element ID, IP ID, face ID that point to me]
2020-01-03 00:17:48 +05:30
real(pReal), dimension(:,:), allocatable, protected :: &
2019-06-06 14:38:58 +05:30
geometry_plastic_nonlocal_IPvolume0 !< volume associated with IP (initially!)
2020-01-03 00:17:48 +05:30
real(pReal), dimension(:,:,:), allocatable, protected :: &
2019-06-06 14:38:58 +05:30
geometry_plastic_nonlocal_IParea0 !< area of interface to neighboring IP (initially!)
2020-01-03 00:17:48 +05:30
real(pReal), dimension(:,:,:,:), allocatable, protected :: &
2019-06-06 14:38:58 +05:30
geometry_plastic_nonlocal_IPareaNormal0 !< area normal of interface to neighboring IP (initially!)
contains
!---------------------------------------------------------------------------------------------------
!> @brief Set the integration point (IP) neighborhood
!> @details: The IP neighborhood for element ID (last index), IP ID (second but last index) and
! face ID (second index) gives the element ID (1 @ first index), IP ID (2 @ first index)
! and face ID (3 @ first index).
! A triangle (2D) has 3 faces, a quadrilateral (2D) had 4 faces, a tetrahedron (3D) has
! 4 faces, and a hexahedron (3D) has 6 faces.
!---------------------------------------------------------------------------------------------------
2019-06-07 13:50:56 +05:30
subroutine geometry_plastic_nonlocal_setIPneighborhood(IPneighborhood)
2019-05-14 15:22:28 +05:30
integer, dimension(:,:,:,:), intent(in) :: IPneighborhood
geometry_plastic_nonlocal_IPneighborhood = IPneighborhood
geometry_plastic_nonlocal_nIPneighbors = size(IPneighborhood,2)
2019-05-14 15:22:28 +05:30
2019-06-07 13:50:56 +05:30
end subroutine geometry_plastic_nonlocal_setIPneighborhood
2019-05-14 15:22:28 +05:30
!---------------------------------------------------------------------------------------------------
!> @brief Set the initial volume associated with an integration point
!---------------------------------------------------------------------------------------------------
2019-06-07 13:50:56 +05:30
subroutine geometry_plastic_nonlocal_setIPvolume(IPvolume)
2019-05-14 15:22:28 +05:30
real(pReal), dimension(:,:), intent(in) :: IPvolume
2019-06-06 14:38:58 +05:30
geometry_plastic_nonlocal_IPvolume0 = IPvolume
2019-05-14 15:22:28 +05:30
2019-06-07 13:50:56 +05:30
end subroutine geometry_plastic_nonlocal_setIPvolume
2019-05-14 15:22:28 +05:30
!---------------------------------------------------------------------------------------------------
!> @brief Set the initial areas of the unit triangle/unit quadrilateral/tetrahedron/hexahedron
! encompassing an integration point
!---------------------------------------------------------------------------------------------------
2019-06-07 13:50:56 +05:30
subroutine geometry_plastic_nonlocal_setIParea(IParea)
real(pReal), dimension(:,:,:), intent(in) :: IParea
geometry_plastic_nonlocal_IParea0 = IParea
2019-06-07 13:50:56 +05:30
end subroutine geometry_plastic_nonlocal_setIParea
!---------------------------------------------------------------------------------------------------
!> @brief Set the direction normal of the areas of the triangle/quadrilateral/tetrahedron/hexahedron
! encompassing an integration point
!---------------------------------------------------------------------------------------------------
2019-06-07 13:50:56 +05:30
subroutine geometry_plastic_nonlocal_setIPareaNormal(IPareaNormal)
real(pReal), dimension(:,:,:,:), intent(in) :: IPareaNormal
geometry_plastic_nonlocal_IPareaNormal0 = IPareaNormal
2019-06-07 13:50:56 +05:30
end subroutine geometry_plastic_nonlocal_setIPareaNormal
2019-06-07 13:50:56 +05:30
!---------------------------------------------------------------------------------------------------
2020-01-03 00:17:48 +05:30
!> @brief Free memory used by variables only needed by plastic_nonlocal
2019-06-07 13:50:56 +05:30
!---------------------------------------------------------------------------------------------------
subroutine geometry_plastic_nonlocal_disable
if(allocated(geometry_plastic_nonlocal_IPneighborhood)) &
deallocate(geometry_plastic_nonlocal_IPneighborhood)
if(allocated(geometry_plastic_nonlocal_IPvolume0)) &
deallocate(geometry_plastic_nonlocal_IPvolume0)
if(allocated(geometry_plastic_nonlocal_IParea0)) &
deallocate(geometry_plastic_nonlocal_IParea0)
if(allocated(geometry_plastic_nonlocal_IPareaNormal0)) &
deallocate(geometry_plastic_nonlocal_IPareaNormal0)
end subroutine geometry_plastic_nonlocal_disable
!---------------------------------------------------------------------------------------------------
2020-01-03 00:17:48 +05:30
!> @brief Write geometry data to results file
!---------------------------------------------------------------------------------------------------
subroutine geometry_plastic_nonlocal_results
2019-10-17 03:51:48 +05:30
2019-10-17 11:18:57 +05:30
integer, dimension(:), allocatable :: shp
call results_openJobFile
2019-10-17 03:51:48 +05:30
writeVolume: block
real(pReal), dimension(:), allocatable :: temp
2019-10-17 11:18:57 +05:30
shp = shape(geometry_plastic_nonlocal_IPvolume0)
temp = reshape(geometry_plastic_nonlocal_IPvolume0,[shp(1)*shp(2)])
2019-10-17 03:51:48 +05:30
call results_writeDataset('geometry',temp,'v_0',&
'initial cell volume','m³')
end block writeVolume
2019-10-17 11:18:57 +05:30
writeAreas: block
2019-10-17 03:51:48 +05:30
real(pReal), dimension(:,:), allocatable :: temp
2019-10-17 11:18:57 +05:30
shp = shape(geometry_plastic_nonlocal_IParea0)
temp = reshape(geometry_plastic_nonlocal_IParea0,[shp(1),shp(2)*shp(3)])
2019-10-17 03:51:48 +05:30
call results_writeDataset('geometry',temp,'a_0',&
'initial cell face area','m²')
2019-10-17 11:18:57 +05:30
end block writeAreas
writeNormals: block
real(pReal), dimension(:,:,:), allocatable :: temp
shp = shape(geometry_plastic_nonlocal_IPareaNormal0)
temp = reshape(geometry_plastic_nonlocal_IPareaNormal0,[shp(1),shp(2),shp(3)*shp(4)])
call results_writeDataset('geometry',temp,'n_0',&
'initial cell face normals','-',transposed=.false.)
end block writeNormals
2019-10-17 03:51:48 +05:30
call results_closeJobFile
end subroutine geometry_plastic_nonlocal_results
end module geometry_plastic_nonlocal