DAMASK_EICMD/src/grid/grid_thermal_spectral.f90

314 lines
15 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
2019-03-12 16:06:18 +05:30
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
!> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
!> @author Shaokang Zhang, Max-Planck-Institut für Eisenforschung GmbH
!> @brief Spectral solver for thermal conduction
!--------------------------------------------------------------------------------------------------
2019-03-12 04:07:06 +05:30
module grid_thermal_spectral
#include <petsc/finclude/petscsnes.h>
#include <petsc/finclude/petscdmda.h>
2019-03-12 16:06:18 +05:30
use PETScdmda
use PETScsnes
2019-06-10 00:50:38 +05:30
use prec
use spectral_utilities
2019-09-28 03:04:34 +05:30
use mesh_grid
2019-06-10 00:50:38 +05:30
use thermal_conduction
use material
use numerics
2019-03-12 16:06:18 +05:30
implicit none
private
2019-06-10 00:50:38 +05:30
!--------------------------------------------------------------------------------------------------
! derived types
2019-03-12 16:06:18 +05:30
type(tSolutionParams), private :: params
!--------------------------------------------------------------------------------------------------
! PETSc data
2019-03-12 16:06:18 +05:30
SNES, private :: thermal_snes
Vec, private :: solution_vec
PetscInt, private :: xstart, xend, ystart, yend, zstart, zend
real(pReal), private, dimension(:,:,:), allocatable :: &
temperature_current, & !< field of current temperature
temperature_lastInc, & !< field of previous temperature
temperature_stagInc !< field of staggered temperature
!--------------------------------------------------------------------------------------------------
! reference diffusion tensor, mobility etc.
2019-03-12 16:06:18 +05:30
integer, private :: totalIter = 0 !< total iteration in current increment
real(pReal), dimension(3,3), private :: D_ref
real(pReal), private :: mobility_ref
public :: &
grid_thermal_spectral_init, &
grid_thermal_spectral_solution, &
grid_thermal_spectral_forward
private :: &
formResidual
contains
!--------------------------------------------------------------------------------------------------
2019-03-12 10:23:12 +05:30
!> @brief allocates all neccessary fields and fills them with data
! ToDo: Restart not implemented
!--------------------------------------------------------------------------------------------------
2019-03-12 04:07:06 +05:30
subroutine grid_thermal_spectral_init
2019-03-09 14:24:33 +05:30
2019-03-12 16:06:18 +05:30
PetscInt, dimension(worldsize) :: localK
2019-03-09 14:24:33 +05:30
integer :: i, j, k, cell
DM :: thermal_grid
2019-03-12 04:07:06 +05:30
PetscScalar, dimension(:,:,:), pointer :: x_scal
2019-03-09 14:24:33 +05:30
PetscErrorCode :: ierr
2019-03-12 04:07:06 +05:30
write(6,'(/,a)') ' <<<+- grid_thermal_spectral init -+>>>'
2019-03-09 15:32:12 +05:30
write(6,'(/,a)') ' Shanthraj et al., Handbook of Mechanics of Materials, 2019'
write(6,'(a)') ' https://doi.org/10.1007/978-981-10-6855-3_80'
2019-03-12 04:07:06 +05:30
!--------------------------------------------------------------------------------------------------
! set default and user defined options for PETSc
call PETScOptionsInsertString(PETSC_NULL_OPTIONS,'-thermal_snes_type ngmres',ierr)
CHKERRQ(ierr)
call PETScOptionsInsertString(PETSC_NULL_OPTIONS,trim(petsc_options),ierr)
CHKERRQ(ierr)
!--------------------------------------------------------------------------------------------------
! initialize solver specific parts of PETSc
2019-03-09 14:24:33 +05:30
call SNESCreate(PETSC_COMM_WORLD,thermal_snes,ierr); CHKERRQ(ierr)
call SNESSetOptionsPrefix(thermal_snes,'thermal_',ierr);CHKERRQ(ierr)
localK = 0
localK(worldrank+1) = grid3
call MPI_Allreduce(MPI_IN_PLACE,localK,worldsize,MPI_INTEGER,MPI_SUM,PETSC_COMM_WORLD,ierr)
call DMDACreate3D(PETSC_COMM_WORLD, &
DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, & ! cut off stencil at boundary
DMDA_STENCIL_BOX, & ! Moore (26) neighborhood around central point
grid(1),grid(2),grid(3), & ! global grid
1, 1, worldsize, &
2019-03-12 16:06:18 +05:30
1, 0, & ! #dof (thermal phase field), ghost boundary width (domain overlap)
[grid(1)],[grid(2)],localK, & ! local grid
thermal_grid,ierr) ! handle, error
2019-03-09 14:24:33 +05:30
CHKERRQ(ierr)
call SNESSetDM(thermal_snes,thermal_grid,ierr); CHKERRQ(ierr) ! connect snes to da
call DMsetFromOptions(thermal_grid,ierr); CHKERRQ(ierr)
call DMsetUp(thermal_grid,ierr); CHKERRQ(ierr)
2019-03-12 04:07:06 +05:30
call DMCreateGlobalVector(thermal_grid,solution_vec,ierr); CHKERRQ(ierr) ! global solution vector (grid x 1, i.e. every def grad tensor)
2019-03-12 10:23:12 +05:30
call DMDASNESSetFunctionLocal(thermal_grid,INSERT_VALUES,formResidual,PETSC_NULL_SNES,ierr) ! residual vector of same shape as solution vector
2019-03-09 14:24:33 +05:30
CHKERRQ(ierr)
call SNESSetFromOptions(thermal_snes,ierr); CHKERRQ(ierr) ! pull it all together with additional CLI arguments
!--------------------------------------------------------------------------------------------------
! init fields
2019-03-09 14:24:33 +05:30
call DMDAGetCorners(thermal_grid,xstart,ystart,zstart,xend,yend,zend,ierr)
CHKERRQ(ierr)
xend = xstart + xend - 1
yend = ystart + yend - 1
zend = zstart + zend - 1
allocate(temperature_current(grid(1),grid(2),grid3), source=0.0_pReal)
allocate(temperature_lastInc(grid(1),grid(2),grid3), source=0.0_pReal)
allocate(temperature_stagInc(grid(1),grid(2),grid3), source=0.0_pReal)
cell = 0
2019-03-12 16:06:18 +05:30
do k = 1, grid3; do j = 1, grid(2); do i = 1,grid(1)
2019-03-09 14:24:33 +05:30
cell = cell + 1
2019-03-10 15:53:39 +05:30
temperature_current(i,j,k) = temperature(material_homogenizationAt(cell))% &
p(thermalMapping(material_homogenizationAt(cell))%p(1,cell))
2019-03-09 14:24:33 +05:30
temperature_lastInc(i,j,k) = temperature_current(i,j,k)
temperature_stagInc(i,j,k) = temperature_current(i,j,k)
enddo; enddo; enddo
2019-03-12 04:07:06 +05:30
call DMDAVecGetArrayF90(thermal_grid,solution_vec,x_scal,ierr); CHKERRQ(ierr) !< get the data out of PETSc to work with
2019-03-09 14:24:33 +05:30
x_scal(xstart:xend,ystart:yend,zstart:zend) = temperature_current
2019-03-12 04:07:06 +05:30
call DMDAVecRestoreArrayF90(thermal_grid,solution_vec,x_scal,ierr); CHKERRQ(ierr)
2016-06-27 21:20:43 +05:30
!--------------------------------------------------------------------------------------------------
! thermal reference diffusion update
2019-03-09 14:24:33 +05:30
cell = 0
D_ref = 0.0_pReal
mobility_ref = 0.0_pReal
do k = 1, grid3; do j = 1, grid(2); do i = 1,grid(1)
cell = cell + 1
D_ref = D_ref + thermal_conduction_getConductivity33(1,cell)
mobility_ref = mobility_ref + thermal_conduction_getMassDensity(1,cell)* &
thermal_conduction_getSpecificHeat(1,cell)
enddo; enddo; enddo
D_ref = D_ref*wgt
call MPI_Allreduce(MPI_IN_PLACE,D_ref,9,MPI_DOUBLE,MPI_SUM,PETSC_COMM_WORLD,ierr)
mobility_ref = mobility_ref*wgt
call MPI_Allreduce(MPI_IN_PLACE,mobility_ref,1,MPI_DOUBLE,MPI_SUM,PETSC_COMM_WORLD,ierr)
2019-03-12 04:07:06 +05:30
end subroutine grid_thermal_spectral_init
2019-03-12 16:06:18 +05:30
!--------------------------------------------------------------------------------------------------
2016-06-27 21:20:43 +05:30
!> @brief solution for the spectral thermal scheme with internal iterations
!--------------------------------------------------------------------------------------------------
2019-09-23 19:20:25 +05:30
function grid_thermal_spectral_solution(timeinc,timeinc_old) result(solution)
2019-03-09 14:24:33 +05:30
real(pReal), intent(in) :: &
2019-03-12 04:07:06 +05:30
timeinc, & !< increment in time for current solution
2019-09-23 19:20:25 +05:30
timeinc_old !< increment in time of last increment
2019-03-09 14:24:33 +05:30
integer :: i, j, k, cell
2019-03-12 04:07:06 +05:30
type(tSolutionState) :: solution
2019-03-09 14:24:33 +05:30
PetscInt :: position
PetscReal :: minTemperature, maxTemperature, stagNorm, solnNorm
2016-06-27 21:20:43 +05:30
2019-03-09 14:24:33 +05:30
PetscErrorCode :: ierr
SNESConvergedReason :: reason
2016-06-27 21:20:43 +05:30
2019-03-12 04:07:06 +05:30
solution%converged =.false.
!--------------------------------------------------------------------------------------------------
! set module wide availabe data
2019-03-09 14:24:33 +05:30
params%timeinc = timeinc
params%timeincOld = timeinc_old
2019-03-12 04:07:06 +05:30
call SNESSolve(thermal_snes,PETSC_NULL_VEC,solution_vec,ierr); CHKERRQ(ierr)
2019-03-09 14:24:33 +05:30
call SNESGetConvergedReason(thermal_snes,reason,ierr); CHKERRQ(ierr)
2019-03-09 14:24:33 +05:30
if (reason < 1) then
2019-03-12 04:07:06 +05:30
solution%converged = .false.
solution%iterationsNeeded = itmax
2019-03-09 14:24:33 +05:30
else
2019-03-12 04:07:06 +05:30
solution%converged = .true.
solution%iterationsNeeded = totalIter
2019-03-09 14:24:33 +05:30
endif
stagNorm = maxval(abs(temperature_current - temperature_stagInc))
solnNorm = maxval(abs(temperature_current))
call MPI_Allreduce(MPI_IN_PLACE,stagNorm,1,MPI_DOUBLE,MPI_MAX,PETSC_COMM_WORLD,ierr)
call MPI_Allreduce(MPI_IN_PLACE,solnNorm,1,MPI_DOUBLE,MPI_MAX,PETSC_COMM_WORLD,ierr)
temperature_stagInc = temperature_current
2019-04-04 03:12:30 +05:30
solution%stagConverged = stagNorm < max(err_thermal_tolAbs, err_thermal_tolRel*solnNorm)
!--------------------------------------------------------------------------------------------------
! updating thermal state
2019-03-09 14:24:33 +05:30
cell = 0
do k = 1, grid3; do j = 1, grid(2); do i = 1,grid(1)
cell = cell + 1
call thermal_conduction_putTemperatureAndItsRate(temperature_current(i,j,k), &
(temperature_current(i,j,k)-temperature_lastInc(i,j,k))/params%timeinc, &
1,cell)
enddo; enddo; enddo
2019-03-12 04:07:06 +05:30
call VecMin(solution_vec,position,minTemperature,ierr); CHKERRQ(ierr)
call VecMax(solution_vec,position,maxTemperature,ierr); CHKERRQ(ierr)
if (solution%converged) &
2019-03-09 14:24:33 +05:30
write(6,'(/,a)') ' ... thermal conduction converged ..................................'
write(6,'(/,a,f8.4,2x,f8.4,2x,f8.4,/)',advance='no') ' Minimum|Maximum|Delta Temperature / K = ',&
minTemperature, maxTemperature, stagNorm
write(6,'(/,a)') ' ==========================================================================='
flush(6)
2019-03-12 04:07:06 +05:30
end function grid_thermal_spectral_solution
2019-03-12 10:23:12 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief forwarding routine
!--------------------------------------------------------------------------------------------------
subroutine grid_thermal_spectral_forward(cutBack)
2019-06-10 00:50:38 +05:30
logical, intent(in) :: cutBack
2019-03-12 16:06:18 +05:30
integer :: i, j, k, cell
DM :: dm_local
PetscScalar, dimension(:,:,:), pointer :: x_scal
PetscErrorCode :: ierr
if (cutBack) then
temperature_current = temperature_lastInc
temperature_stagInc = temperature_lastInc
2019-03-12 10:23:12 +05:30
!--------------------------------------------------------------------------------------------------
! reverting thermal field state
2019-03-12 16:06:18 +05:30
cell = 0
call SNESGetDM(thermal_snes,dm_local,ierr); CHKERRQ(ierr)
call DMDAVecGetArrayF90(dm_local,solution_vec,x_scal,ierr); CHKERRQ(ierr) !< get the data out of PETSc to work with
x_scal(xstart:xend,ystart:yend,zstart:zend) = temperature_current
call DMDAVecRestoreArrayF90(dm_local,solution_vec,x_scal,ierr); CHKERRQ(ierr)
do k = 1, grid3; do j = 1, grid(2); do i = 1,grid(1)
cell = cell + 1
call thermal_conduction_putTemperatureAndItsRate(temperature_current(i,j,k), &
(temperature_current(i,j,k) - &
temperature_lastInc(i,j,k))/params%timeinc, &
1,cell)
enddo; enddo; enddo
else
2019-03-12 10:23:12 +05:30
!--------------------------------------------------------------------------------------------------
! update rate and forward last inc
2019-03-12 16:06:18 +05:30
temperature_lastInc = temperature_current
cell = 0
D_ref = 0.0_pReal
mobility_ref = 0.0_pReal
do k = 1, grid3; do j = 1, grid(2); do i = 1,grid(1)
cell = cell + 1
D_ref = D_ref + thermal_conduction_getConductivity33(1,cell)
mobility_ref = mobility_ref + thermal_conduction_getMassDensity(1,cell)* &
thermal_conduction_getSpecificHeat(1,cell)
enddo; enddo; enddo
D_ref = D_ref*wgt
call MPI_Allreduce(MPI_IN_PLACE,D_ref,9,MPI_DOUBLE,MPI_SUM,PETSC_COMM_WORLD,ierr)
mobility_ref = mobility_ref*wgt
call MPI_Allreduce(MPI_IN_PLACE,mobility_ref,1,MPI_DOUBLE,MPI_SUM,PETSC_COMM_WORLD,ierr)
endif
2019-03-12 10:23:12 +05:30
end subroutine grid_thermal_spectral_forward
!--------------------------------------------------------------------------------------------------
!> @brief forms the spectral thermal residual vector
!--------------------------------------------------------------------------------------------------
2019-03-12 10:23:12 +05:30
subroutine formResidual(in,x_scal,f_scal,dummy,ierr)
2019-03-12 16:06:18 +05:30
DMDALocalInfo, dimension(DMDA_LOCAL_INFO_SIZE) :: &
in
PetscScalar, dimension( &
XG_RANGE,YG_RANGE,ZG_RANGE), intent(in) :: &
x_scal
PetscScalar, dimension( &
X_RANGE,Y_RANGE,Z_RANGE), intent(out) :: &
f_scal
PetscObject :: dummy
PetscErrorCode :: ierr
integer :: i, j, k, cell
real(pReal) :: Tdot, dTdot_dT
2019-03-12 16:06:18 +05:30
temperature_current = x_scal
!--------------------------------------------------------------------------------------------------
! evaluate polarization field
2019-03-12 16:06:18 +05:30
scalarField_real = 0.0_pReal
scalarField_real(1:grid(1),1:grid(2),1:grid3) = temperature_current
call utilities_FFTscalarForward
call utilities_fourierScalarGradient !< calculate gradient of damage field
call utilities_FFTvectorBackward
cell = 0
do k = 1, grid3; do j = 1, grid(2); do i = 1,grid(1)
cell = cell + 1
vectorField_real(1:3,i,j,k) = matmul(thermal_conduction_getConductivity33(1,cell) - D_ref, &
2019-03-12 16:06:18 +05:30
vectorField_real(1:3,i,j,k))
enddo; enddo; enddo
call utilities_FFTvectorForward
call utilities_fourierVectorDivergence !< calculate damage divergence in fourier field
call utilities_FFTscalarBackward
cell = 0
do k = 1, grid3; do j = 1, grid(2); do i = 1,grid(1)
cell = cell + 1
call thermal_conduction_getSourceAndItsTangent(Tdot, dTdot_dT, temperature_current(i,j,k), 1, cell)
scalarField_real(i,j,k) = params%timeinc*scalarField_real(i,j,k) + &
params%timeinc*Tdot + &
thermal_conduction_getMassDensity (1,cell)* &
thermal_conduction_getSpecificHeat(1,cell)*(temperature_lastInc(i,j,k) - &
temperature_current(i,j,k)) + &
mobility_ref*temperature_current(i,j,k)
enddo; enddo; enddo
!--------------------------------------------------------------------------------------------------
! convolution of damage field with green operator
2019-03-12 16:06:18 +05:30
call utilities_FFTscalarForward
call utilities_fourierGreenConvolution(D_ref, mobility_ref, params%timeinc)
call utilities_FFTscalarBackward
!--------------------------------------------------------------------------------------------------
! constructing residual
2019-03-12 16:06:18 +05:30
f_scal = temperature_current - scalarField_real(1:grid(1),1:grid(2),1:grid3)
2019-03-12 10:23:12 +05:30
end subroutine formResidual
2019-03-12 04:07:06 +05:30
end module grid_thermal_spectral