2016-01-17 20:33:54 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @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
|
2019-06-06 17:29:16 +05:30
|
|
|
use prec
|
|
|
|
use numerics
|
|
|
|
use debug
|
|
|
|
use config
|
|
|
|
use FEsolving
|
|
|
|
use math
|
2019-09-22 19:23:03 +05:30
|
|
|
use rotations
|
2019-06-06 17:29:16 +05:30
|
|
|
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
|
2019-06-06 17:29:16 +05:30
|
|
|
|
|
|
|
implicit none
|
|
|
|
private
|
|
|
|
|
|
|
|
public :: &
|
2019-10-28 17:59:32 +05:30
|
|
|
CPFEM_forward, &
|
2019-06-06 17:29:16 +05:30
|
|
|
CPFEM_initAll, &
|
2019-10-24 09:46:42 +05:30
|
|
|
CPFEM_results, &
|
|
|
|
CPFEM_restartWrite
|
2019-05-15 02:42:32 +05:30
|
|
|
|
2016-01-17 20:33:54 +05:30
|
|
|
contains
|
|
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-12-07 19:50:04 +05:30
|
|
|
!> @brief call all module initializations
|
2016-01-17 20:33:54 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-06-11 19:46:10 +05:30
|
|
|
subroutine CPFEM_initAll
|
2016-01-17 20:33:54 +05:30
|
|
|
|
2019-06-15 19:18:47 +05:30
|
|
|
call DAMASK_interface_init ! Spectral and FEM interface to commandline
|
2019-06-11 19:46:10 +05:30
|
|
|
call prec_init
|
|
|
|
call IO_init
|
2016-01-17 20:33:54 +05:30
|
|
|
#ifdef FEM
|
2019-06-11 19:46:10 +05:30
|
|
|
call FEM_Zoo_init
|
2016-01-17 20:33:54 +05:30
|
|
|
#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
|
2016-01-17 20:33:54 +05:30
|
|
|
|
|
|
|
end subroutine CPFEM_initAll
|
|
|
|
|
2019-12-07 19:50:04 +05:30
|
|
|
|
2016-01-17 20:33:54 +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
|
|
|
|
2019-06-15 19:18:47 +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'
|
2019-12-09 01:22:05 +05:30
|
|
|
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'
|
2019-12-09 01:22:05 +05:30
|
|
|
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
|
2016-01-17 20:33:54 +05:30
|
|
|
|
|
|
|
end subroutine CPFEM_init
|
|
|
|
|
2018-12-05 04:25:39 +05:30
|
|
|
|
2016-01-17 20:33:54 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-10-24 09:46:42 +05:30
|
|
|
!> @brief Forward data after successful increment.
|
2019-10-28 18:08:30 +05:30
|
|
|
! ToDo: Any guessing for the current states possible?
|
2016-01-17 20:33:54 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-10-28 17:59:32 +05:30
|
|
|
subroutine CPFEM_forward
|
2019-10-24 02:36:47 +05:30
|
|
|
|
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
|
2019-06-15 19:18:47 +05:30
|
|
|
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-10-24 09:46:42 +05:30
|
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-12-07 19:50:04 +05:30
|
|
|
!> @brief Write current restart information (Field and constitutive data) to file.
|
2019-10-24 09:46:42 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
subroutine CPFEM_restartWrite
|
|
|
|
|
2019-12-07 19:50:04 +05:30
|
|
|
integer :: i
|
|
|
|
integer(HID_T) :: fileHandle, groupHandle
|
|
|
|
character(len=pStringLen) :: fileName, datasetName
|
2019-10-24 16:36:42 +05:30
|
|
|
|
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)
|
2019-10-24 09:46:42 +05:30
|
|
|
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)
|
2019-10-24 09:46:42 +05:30
|
|
|
enddo
|
2019-12-07 19:50:04 +05:30
|
|
|
call HDF5_closeGroup(groupHandle)
|
2019-06-11 19:46:10 +05:30
|
|
|
|
2019-10-24 09:46:42 +05:30
|
|
|
call HDF5_closeFile(fileHandle)
|
2019-06-11 19:46:10 +05:30
|
|
|
|
2019-10-24 09:46:42 +05:30
|
|
|
end subroutine CPFEM_restartWrite
|
2016-01-17 20:33:54 +05:30
|
|
|
|
2019-05-05 15:36:55 +05:30
|
|
|
|
2018-12-05 04:25:39 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-10-24 09:46:42 +05:30
|
|
|
!> @brief Trigger writing of results.
|
2018-12-05 04:25:39 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-12-17 20:45:16 +05:30
|
|
|
subroutine CPFEM_results(inc,time)
|
2019-05-15 02:42:32 +05:30
|
|
|
|
2019-06-15 19:18:47 +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
|
2020-01-10 06:15:00 +05:30
|
|
|
call results_finalizeIncrement
|
2019-06-11 19:46:10 +05:30
|
|
|
call results_closeJobFile
|
2018-12-05 04:25:39 +05:30
|
|
|
|
|
|
|
end subroutine CPFEM_results
|
|
|
|
|
2016-01-17 20:33:54 +05:30
|
|
|
end module CPFEM2
|