DAMASK_EICMD/src/CPFEM2.f90

210 lines
6.6 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
!> @brief needs a good name and description
!--------------------------------------------------------------------------------------------------
module CPFEM2
use prec
use numerics
use debug
use config
use FEsolving
use math
2019-09-22 19:23:03 +05:30
use rotations
use material
use lattice
use IO
use DAMASK_interface
use results
use discretization
use HDF5_utilities
use homogenization
use constitutive
use crystallite
2019-06-11 19:46:10 +05:30
#ifdef FEM
use FEM_Zoo
2019-09-28 03:04:34 +05:30
use mesh
#else
use mesh_grid
2019-06-11 19:46:10 +05:30
#endif
implicit none
private
public :: &
2019-10-28 17:59:32 +05:30
CPFEM_forward, &
CPFEM_initAll, &
CPFEM_results, &
CPFEM_restartWrite
2019-05-15 02:42:32 +05:30
contains
!--------------------------------------------------------------------------------------------------
2019-12-07 19:50:04 +05:30
!> @brief call all module initializations
!--------------------------------------------------------------------------------------------------
2019-06-11 19:46:10 +05:30
subroutine CPFEM_initAll
call DAMASK_interface_init ! Spectral and FEM interface to commandline
2019-06-11 19:46:10 +05:30
call prec_init
call IO_init
#ifdef FEM
2019-06-11 19:46:10 +05:30
call FEM_Zoo_init
#endif
2019-06-11 19:46:10 +05:30
call numerics_init
call debug_init
call config_init
call math_init
2019-09-22 19:23:03 +05:30
call rotations_init
2019-06-11 19:46:10 +05:30
call lattice_init
call HDF5_utilities_init
call results_init
2019-11-24 12:14:17 +05:30
call mesh_init
2019-06-11 19:46:10 +05:30
call material_init
call constitutive_init
call crystallite_init
call homogenization_init
call CPFEM_init
end subroutine CPFEM_initAll
2019-12-07 19:50:04 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief allocate the arrays defined in module CPFEM and initialize them
!--------------------------------------------------------------------------------------------------
subroutine CPFEM_init
2019-05-15 02:42:32 +05:30
2019-12-07 19:50:04 +05:30
integer :: i
integer(HID_T) :: fileHandle, groupHandle
character(len=pStringLen) :: fileName, datasetName
2019-06-11 19:46:10 +05:30
2019-12-07 19:50:04 +05:30
write(6,'(/,a)') ' <<<+- CPFEM init -+>>>'; flush(6)
2019-06-11 19:46:10 +05:30
if (interface_restartInc > 0) then
2019-12-07 19:50:04 +05:30
write(6,'(/,a,i0,a)') ' reading restart information of increment ', interface_restartInc, ' from file'
2019-06-11 19:46:10 +05:30
2019-12-07 19:50:04 +05:30
write(fileName,'(a,i0,a)') trim(getSolverJobName())//'_',worldrank,'.hdf5'
fileHandle = HDF5_openFile(fileName)
2019-10-28 17:47:05 +05:30
call HDF5_read(fileHandle,crystallite_F0, 'F')
call HDF5_read(fileHandle,crystallite_Fp0,'Fp')
call HDF5_read(fileHandle,crystallite_Fi0,'Fi')
call HDF5_read(fileHandle,crystallite_Lp0,'Lp')
call HDF5_read(fileHandle,crystallite_Li0,'Li')
call HDF5_read(fileHandle,crystallite_S0, 'S')
2019-06-11 19:46:10 +05:30
2019-12-07 19:50:04 +05:30
groupHandle = HDF5_openGroup(fileHandle,'constituent')
do i = 1,size(phase_plasticity)
write(datasetName,'(i0,a)') i,'_omega_plastic'
call HDF5_read(groupHandle,plasticState(i)%state0,datasetName)
2019-06-11 19:46:10 +05:30
enddo
2019-12-07 19:50:04 +05:30
call HDF5_closeGroup(groupHandle)
groupHandle = HDF5_openGroup(fileHandle,'materialpoint')
do i = 1, material_Nhomogenization
write(datasetName,'(i0,a)') i,'_omega_homogenization'
call HDF5_read(groupHandle,homogState(i)%state0,datasetName)
2019-06-11 19:46:10 +05:30
enddo
2019-12-07 19:50:04 +05:30
call HDF5_closeGroup(groupHandle)
2019-06-11 19:46:10 +05:30
call HDF5_closeFile(fileHandle)
endif
end subroutine CPFEM_init
2018-12-05 04:25:39 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Forward data after successful increment.
! ToDo: Any guessing for the current states possible?
!--------------------------------------------------------------------------------------------------
2019-10-28 17:59:32 +05:30
subroutine CPFEM_forward
2019-12-07 19:50:04 +05:30
integer :: i, j
2019-06-11 19:46:10 +05:30
if (iand(debug_level(debug_CPFEM), debug_levelBasic) /= 0) &
write(6,'(a)') '<< CPFEM >> aging states'
crystallite_F0 = crystallite_partionedF
crystallite_Fp0 = crystallite_Fp
crystallite_Lp0 = crystallite_Lp
crystallite_Fi0 = crystallite_Fi
crystallite_Li0 = crystallite_Li
crystallite_S0 = crystallite_S
2019-06-11 19:46:10 +05:30
do i = 1, size(plasticState)
plasticState(i)%state0 = plasticState(i)%state
enddo
do i = 1, size(sourceState)
2019-12-07 19:50:04 +05:30
do j = 1,phase_Nsources(i)
sourceState(i)%p(j)%state0 = sourceState(i)%p(j)%state
2019-06-11 19:46:10 +05:30
enddo; enddo
2019-12-07 19:50:04 +05:30
do i = 1, material_Nhomogenization
homogState (i)%state0 = homogState (i)%state
thermalState(i)%state0 = thermalState(i)%state
damageState (i)%state0 = damageState (i)%state
2019-06-11 19:46:10 +05:30
enddo
2019-10-28 17:59:32 +05:30
end subroutine CPFEM_forward
!--------------------------------------------------------------------------------------------------
2019-12-07 19:50:04 +05:30
!> @brief Write current restart information (Field and constitutive data) to file.
!--------------------------------------------------------------------------------------------------
subroutine CPFEM_restartWrite
2019-12-07 19:50:04 +05:30
integer :: i
integer(HID_T) :: fileHandle, groupHandle
character(len=pStringLen) :: fileName, datasetName
2019-12-07 19:50:04 +05:30
write(6,'(a)') ' writing field and constitutive data required for restart to file';flush(6)
2019-06-11 19:46:10 +05:30
2019-12-07 19:50:04 +05:30
write(fileName,'(a,i0,a)') trim(getSolverJobName())//'_',worldrank,'.hdf5'
fileHandle = HDF5_openFile(fileName,'a')
2019-06-11 19:46:10 +05:30
2019-10-28 17:47:05 +05:30
call HDF5_write(fileHandle,crystallite_partionedF,'F')
call HDF5_write(fileHandle,crystallite_Fp, 'Fp')
call HDF5_write(fileHandle,crystallite_Fi, 'Fi')
call HDF5_write(fileHandle,crystallite_Lp, 'Lp')
call HDF5_write(fileHandle,crystallite_Li, 'Li')
call HDF5_write(fileHandle,crystallite_S, 'S')
2019-06-11 19:46:10 +05:30
2019-12-07 19:50:04 +05:30
groupHandle = HDF5_addGroup(fileHandle,'constituent')
do i = 1,size(phase_plasticity)
write(datasetName,'(i0,a)') i,'_omega_plastic'
call HDF5_write(groupHandle,plasticState(i)%state,datasetName)
enddo
2019-12-07 19:50:04 +05:30
call HDF5_closeGroup(groupHandle)
2019-06-11 19:46:10 +05:30
2019-12-07 19:50:04 +05:30
groupHandle = HDF5_addGroup(fileHandle,'materialpoint')
do i = 1, material_Nhomogenization
write(datasetName,'(i0,a)') i,'_omega_homogenization'
call HDF5_write(groupHandle,homogState(i)%state,datasetName)
enddo
2019-12-07 19:50:04 +05:30
call HDF5_closeGroup(groupHandle)
2019-06-11 19:46:10 +05:30
call HDF5_closeFile(fileHandle)
2019-06-11 19:46:10 +05:30
end subroutine CPFEM_restartWrite
2018-12-05 04:25:39 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Trigger writing of results.
2018-12-05 04:25:39 +05:30
!--------------------------------------------------------------------------------------------------
subroutine CPFEM_results(inc,time)
2019-05-15 02:42:32 +05:30
integer, intent(in) :: inc
real(pReal), intent(in) :: time
2019-06-11 19:46:10 +05:30
call results_openJobFile
call results_addIncrement(inc,time)
call constitutive_results
call crystallite_results
call homogenization_results
call discretization_results
2019-12-07 19:50:04 +05:30
call results_removeLink('current') ! ToDo: put this into closeJobFile?
2019-06-11 19:46:10 +05:30
call results_closeJobFile
2018-12-05 04:25:39 +05:30
end subroutine CPFEM_results
end module CPFEM2