DAMASK_EICMD/src/discretization.f90

121 lines
4.2 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @brief spatial discretization
2019-06-15 19:10:43 +05:30
!> @details serves as an abstraction layer between the different solvers and DAMASK
!--------------------------------------------------------------------------------------------------
module discretization
use prec
use results
implicit none(type,external)
private
2021-05-22 20:51:07 +05:30
2019-06-07 02:19:17 +05:30
integer, public, protected :: &
discretization_nIPs, &
discretization_Nelems, &
discretization_Ncells
2021-05-22 20:51:07 +05:30
2019-06-07 02:19:17 +05:30
integer, public, protected, dimension(:), allocatable :: &
2021-05-22 20:51:07 +05:30
discretization_materialAt !ToDo: discretization_ID_material
2021-05-22 20:51:07 +05:30
real(pReal), public, protected, dimension(:,:), allocatable :: &
discretization_IPcoords0, &
discretization_IPcoords, &
2020-01-27 01:45:21 +05:30
discretization_NodeCoords0, &
discretization_NodeCoords
2021-05-22 20:51:07 +05:30
integer :: &
2020-03-10 02:50:33 +05:30
discretization_sharedNodesBegin
public :: &
discretization_init, &
discretization_results, &
2019-09-28 02:37:34 +05:30
discretization_setIPcoords, &
discretization_setNodeCoords
contains
2021-05-22 20:51:07 +05:30
2019-06-15 19:10:43 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief stores the relevant information in globally accesible variables
!--------------------------------------------------------------------------------------------------
2020-10-01 16:13:05 +05:30
subroutine discretization_init(materialAt,&
IPcoords0,NodeCoords0,&
2020-03-10 02:50:33 +05:30
sharedNodesBegin)
2019-06-15 19:10:43 +05:30
integer, dimension(:), intent(in) :: &
2020-10-01 16:13:05 +05:30
materialAt
real(pReal), dimension(:,:), intent(in) :: &
IPcoords0, &
NodeCoords0
integer, optional, intent(in) :: &
2020-06-21 02:21:00 +05:30
sharedNodesBegin !< index of first node shared among different processes (MPI)
print'(/,1x,a)', '<<<+- discretization init -+>>>'; flush(6)
discretization_Nelems = size(materialAt,1)
discretization_nIPs = size(IPcoords0,2)/discretization_Nelems
discretization_Ncells = discretization_Nelems*discretization_nIPs
2021-05-22 20:51:07 +05:30
discretization_materialAt = materialAt
discretization_IPcoords0 = IPcoords0
discretization_IPcoords = IPcoords0
discretization_NodeCoords0 = NodeCoords0
discretization_NodeCoords = NodeCoords0
2021-05-22 20:51:07 +05:30
2022-12-07 22:59:03 +05:30
if (present(sharedNodesBegin)) then
2020-03-10 02:50:33 +05:30
discretization_sharedNodesBegin = sharedNodesBegin
else
2020-03-10 02:50:33 +05:30
discretization_sharedNodesBegin = size(discretization_NodeCoords0,2)
2022-06-09 02:36:01 +05:30
end if
2021-05-22 20:51:07 +05:30
end subroutine discretization_init
2019-06-15 19:10:43 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief write the displacements
!--------------------------------------------------------------------------------------------------
subroutine discretization_results
2019-12-19 00:35:51 +05:30
real(pReal), dimension(:,:), allocatable :: u
2021-05-22 20:51:07 +05:30
call results_closeGroup(results_addGroup('current/geometry'))
2021-05-22 20:51:07 +05:30
u = discretization_NodeCoords (:,:discretization_sharedNodesBegin) &
- discretization_NodeCoords0(:,:discretization_sharedNodesBegin)
2021-06-01 20:35:13 +05:30
call results_writeDataset(u,'current/geometry','u_n','displacements of the nodes','m')
2021-05-22 20:51:07 +05:30
u = discretization_IPcoords &
- discretization_IPcoords0
2021-06-01 20:35:13 +05:30
call results_writeDataset(u,'current/geometry','u_p','displacements of the materialpoints (cell centers)','m')
2019-12-19 00:35:51 +05:30
end subroutine discretization_results
2019-06-15 19:10:43 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief stores current IP coordinates
!--------------------------------------------------------------------------------------------------
subroutine discretization_setIPcoords(IPcoords)
real(pReal), dimension(:,:), intent(in) :: IPcoords
2021-05-22 20:51:07 +05:30
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
2021-05-22 20:51:07 +05:30
2019-09-28 02:37:34 +05:30
discretization_NodeCoords = NodeCoords
end subroutine discretization_setNodeCoords
end module discretization