2019-06-06 17:29:16 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief spatial discretization
|
2019-06-15 19:10:43 +05:30
|
|
|
!> @details serves as an abstraction layer between the different solvers and DAMASK
|
2019-06-06 17:29:16 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
module discretization
|
|
|
|
|
|
|
|
use prec
|
|
|
|
use results
|
2019-09-20 01:18:04 +05:30
|
|
|
#if defined(PETSc) || defined(DAMASK_HDF5)
|
|
|
|
use HDF5_utilities
|
|
|
|
#endif
|
2019-06-06 17:29:16 +05:30
|
|
|
|
|
|
|
implicit none
|
|
|
|
private
|
|
|
|
|
2019-06-07 02:19:17 +05:30
|
|
|
integer, public, protected :: &
|
2019-06-06 21:58:10 +05:30
|
|
|
discretization_nIP, &
|
|
|
|
discretization_nElem
|
2019-06-07 00:24:19 +05:30
|
|
|
|
2019-06-07 02:19:17 +05:30
|
|
|
integer, public, protected, dimension(:), allocatable :: &
|
2019-06-07 00:24:19 +05:30
|
|
|
discretization_homogenizationAt, &
|
|
|
|
discretization_microstructureAt
|
2019-06-06 21:58:10 +05:30
|
|
|
|
2019-06-07 02:19:17 +05:30
|
|
|
real(pReal), public, protected, dimension(:,:), allocatable :: &
|
2019-06-06 21:58:10 +05:30
|
|
|
discretization_IPcoords0, &
|
|
|
|
discretization_NodeCoords0, &
|
|
|
|
discretization_IPcoords, &
|
|
|
|
discretization_NodeCoords
|
2019-06-06 17:29:16 +05:30
|
|
|
|
|
|
|
public :: &
|
|
|
|
discretization_init, &
|
2019-06-06 21:58:10 +05:30
|
|
|
discretization_results, &
|
2019-09-28 02:37:34 +05:30
|
|
|
discretization_setIPcoords, &
|
|
|
|
discretization_setNodeCoords
|
2019-06-06 17:29:16 +05:30
|
|
|
|
|
|
|
contains
|
|
|
|
|
2019-06-15 19:10:43 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief stores the relevant information in globally accesible variables
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-06-07 00:24:19 +05:30
|
|
|
subroutine discretization_init(homogenizationAt,microstructureAt,IPcoords0,NodeCoords0)
|
2019-06-06 17:29:16 +05:30
|
|
|
|
2019-06-15 19:10:43 +05:30
|
|
|
integer, dimension(:), intent(in) :: &
|
2019-06-07 00:24:19 +05:30
|
|
|
homogenizationAt, &
|
|
|
|
microstructureAt
|
2019-06-06 21:58:10 +05:30
|
|
|
real(pReal), dimension(:,:), intent(in) :: &
|
|
|
|
IPcoords0, &
|
|
|
|
NodeCoords0
|
2019-06-06 17:29:16 +05:30
|
|
|
|
|
|
|
write(6,'(/,a)') ' <<<+- discretization init -+>>>'
|
2019-06-07 00:24:19 +05:30
|
|
|
|
|
|
|
discretization_nElem = size(microstructureAt,1)
|
2019-06-07 02:30:06 +05:30
|
|
|
discretization_nIP = size(IPcoords0,2)/discretization_nElem
|
2019-06-06 17:29:16 +05:30
|
|
|
|
2019-06-07 00:24:19 +05:30
|
|
|
discretization_homogenizationAt = homogenizationAt
|
|
|
|
discretization_microstructureAt = microstructureAt
|
|
|
|
|
2019-06-06 21:58:10 +05:30
|
|
|
discretization_IPcoords0 = IPcoords0
|
|
|
|
discretization_IPcoords = IPcoords0
|
2019-06-07 00:24:19 +05:30
|
|
|
|
2019-06-06 21:58:10 +05:30
|
|
|
discretization_NodeCoords0 = NodeCoords0
|
|
|
|
discretization_NodeCoords = NodeCoords0
|
|
|
|
|
2019-06-06 17:29:16 +05:30
|
|
|
end subroutine discretization_init
|
|
|
|
|
|
|
|
|
2019-06-15 19:10:43 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief write the displacements
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-06-06 17:29:16 +05:30
|
|
|
subroutine discretization_results
|
2019-06-07 00:24:19 +05:30
|
|
|
#if defined(PETSc) || defined(DAMASK_HDF5)
|
2019-06-06 21:58:10 +05:30
|
|
|
real(pReal), dimension(:,:), allocatable :: u
|
|
|
|
|
2019-09-20 01:18:04 +05:30
|
|
|
call HDF5_closeGroup(results_addGroup(trim('current/geometry')))
|
|
|
|
|
|
|
|
u = discretization_NodeCoords - discretization_NodeCoords0
|
|
|
|
call results_writeDataset('current/geometry',u,'u_n','nodal displacements','m')
|
2019-06-06 21:58:10 +05:30
|
|
|
|
2019-09-20 01:18:04 +05:30
|
|
|
u = discretization_IPcoords - discretization_IPcoords0
|
|
|
|
call results_writeDataset('current/geometry',u,'u_c','cell center displacements','m')
|
2019-06-07 00:24:19 +05:30
|
|
|
#endif
|
2019-06-06 17:29:16 +05:30
|
|
|
end subroutine discretization_results
|
|
|
|
|
2019-06-06 21:58:10 +05:30
|
|
|
|
2019-06-15 19:10:43 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief stores current IP coordinates
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-06-06 21:58:10 +05:30
|
|
|
subroutine discretization_setIPcoords(IPcoords)
|
|
|
|
|
|
|
|
real(pReal), dimension(:,:), intent(in) :: IPcoords
|
|
|
|
|
|
|
|
discretization_IPcoords = IPcoords
|
|
|
|
|
|
|
|
end subroutine discretization_setIPcoords
|
|
|
|
|
2019-09-28 02:37:34 +05:30
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief stores current IP coordinates
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
subroutine discretization_setNodeCoords(NodeCoords)
|
|
|
|
|
|
|
|
real(pReal), dimension(:,:), intent(in) :: NodeCoords
|
|
|
|
|
|
|
|
discretization_NodeCoords = NodeCoords
|
|
|
|
|
|
|
|
end subroutine discretization_setNodeCoords
|
|
|
|
|
|
|
|
|
2019-06-06 17:29:16 +05:30
|
|
|
end module discretization
|