added comments and structured the code. Temperature is not longer stored for each point but only to simplify restart behavior as long as it is not fully supported by the constitutive models

This commit is contained in:
Martin Diehl 2012-11-12 14:14:39 +00:00
parent 7974001c9d
commit 70c4e11742
4 changed files with 774 additions and 751 deletions

View File

@ -1,20 +1,18 @@
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! $Id: DAMASK_spectral_SolverAL.f90 1654 2012-08-03 09:25:48Z MPIE\m.diehl $ ! $Id: DAMASK_spectral_solverAL.f90 1654 2012-08-03 09:25:48Z MPIE\m.diehl $
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH !> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH !> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH !> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
!> @brief AL scheme solver !> @brief AL scheme solver
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
module DAMASK_spectral_SolverAL module DAMASK_spectral_solverAL
use prec, only: & use prec, only: &
pInt, & pInt, &
pReal pReal
use math, only: & use math, only: &
math_I3 math_I3
use DAMASK_spectral_utilities, only: &
use DAMASK_spectral_Utilities, only: &
tSolutionState tSolutionState
implicit none implicit none
@ -26,13 +24,14 @@ module DAMASK_spectral_SolverAL
DAMASK_spectral_SolverAL_label = 'al' DAMASK_spectral_SolverAL_label = 'al'
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! derived types ! derived types ToDo: use here the type definition for a full loadcase including mask
type tSolutionParams type tSolutionParams
real(pReal), dimension(3,3) :: P_BC, rotation_BC real(pReal), dimension(3,3) :: P_BC, rotation_BC
real(pReal) :: timeinc real(pReal) :: timeinc
end type tSolutionParams end type tSolutionParams
type(tSolutionParams), private :: params type(tSolutionParams), private :: params
real(pReal), private, dimension(3,3) :: mask_stress = 0.0_pReal
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! PETSc data ! PETSc data
@ -42,27 +41,34 @@ module DAMASK_spectral_SolverAL
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! common pointwise data ! common pointwise data
real(pReal), private, dimension(:,:,:,:,:), allocatable :: F_lastInc, F_lambda_lastInc, F_lambdaDot, Fdot real(pReal), private, dimension(:,:,:,:,:), allocatable :: &
F_lastInc, & !< field of previous compatible deformation gradients
F_lambda_lastInc, & !< field of previous incompatible deformation gradient
Fdot, & !< field of assumed rate of compatible deformation gradient
F_lambdaDot !< field of assumed rate of incopatible deformation gradient
real(pReal), private, dimension(:,:,:,:), allocatable :: coordinates real(pReal), private, dimension(:,:,:,:), allocatable :: coordinates
real(pReal), private, dimension(:,:,:), allocatable :: temperature real(pReal), private :: temperature !< temperature, no spatial quantity at the moment
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! stress, stiffness and compliance average etc. ! stress, stiffness and compliance average etc.
real(pReal), private, dimension(3,3) :: & real(pReal), private, dimension(3,3) :: &
F_aimDot, & F_aimDot, & !< assumed rate of average deformation gradient
F_aim = math_I3, & F_aim = math_I3, & !< current prescribed deformation gradient
F_aim_lastInc = math_I3, & F_aim_lastInc = math_I3, & !< previous average deformation gradient
P_av P_av = 0.0_pReal !< average 1st Piola--Kirchhoff stress
character(len=1024), private :: incInfo character(len=1024), private :: incInfo !< time and increment information
real(pReal), private, dimension(3,3,3,3) :: & real(pReal), private, dimension(3,3,3,3) :: &
C = 0.0_pReal, C_lastInc = 0.0_pReal, & C = 0.0_pReal, & !< current average stiffness
S = 0.0_pReal, & C_lastInc = 0.0_pReal, & !< previous average stiffness
S = 0.0_pReal, & !< current compliance (filled up with zeros)
C_scale = 0.0_pReal, & C_scale = 0.0_pReal, &
S_scale = 0.0_pReal S_scale = 0.0_pReal
real(pReal), private :: err_stress, err_f, err_p real(pReal), private :: &
err_stress, & !< deviation from stress BC
err_f, & !< difference between compatible and incompatible deformation gradient
err_p !< difference of stress resulting from compatible and incompatible F
logical, private :: ForwardData logical, private :: ForwardData
real(pReal), private, dimension(3,3) :: mask_stress = 0.0_pReal
integer(pInt) :: reportIter = 0_pInt integer(pInt) :: reportIter = 0_pInt
contains contains
@ -71,31 +77,24 @@ contains
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
subroutine AL_init() subroutine AL_init()
use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran >4.6 at the moment) use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran >4.6 at the moment)
use IO, only: & use IO, only: &
IO_read_JobBinaryFile, & IO_read_JobBinaryFile, &
IO_write_JobBinaryFile IO_write_JobBinaryFile
use FEsolving, only: & use FEsolving, only: &
restartInc restartInc
use DAMASK_interface, only: & use DAMASK_interface, only: &
getSolverJobName getSolverJobName
use DAMASK_spectral_Utilities, only: & use DAMASK_spectral_Utilities, only: &
Utilities_init, & Utilities_init, &
Utilities_constitutiveResponse, & Utilities_constitutiveResponse, &
Utilities_updateGamma, & Utilities_updateGamma, &
debugRestart debugRestart
use numerics, only: & use numerics, only: &
petsc_options petsc_options
use mesh, only: & use mesh, only: &
res, & res, &
geomdim, & geomdim, &
mesh_NcpElems mesh_NcpElems
use math, only: & use math, only: &
math_invSym3333 math_invSym3333
@ -103,14 +102,13 @@ subroutine AL_init()
#include <finclude/petscdmda.h90> #include <finclude/petscdmda.h90>
#include <finclude/petscsnes.h90> #include <finclude/petscsnes.h90>
integer(pInt) :: i,j,k integer(pInt) :: i,j,k
real(pReal), dimension(:,:,:,:,:), allocatable :: P real(pReal), dimension(3,3, res(1), res(2),res(3)) :: P
PetscErrorCode :: ierr PetscErrorCode :: ierr
PetscObject :: dummy PetscObject :: dummy
PetscScalar, pointer, dimension(:,:,:,:) :: xx_psc, F, F_lambda PetscScalar, pointer, dimension(:,:,:,:) :: xx_psc, F, F_lambda
call Utilities_init() call Utilities_init()
write(6,'(/,a)') ' <<<+- DAMASK_spectral_solverAL init -+>>>' write(6,'(/,a)') ' <<<+- DAMASK_spectral_solverAL init -+>>>'
write(6,'(a)') ' $Id: DAMASK_spectral_SolverAL.f90 1654 2012-08-03 09:25:48Z MPIE\m.diehl $' write(6,'(a)') ' $Id: DAMASK_spectral_SolverAL.f90 1654 2012-08-03 09:25:48Z MPIE\m.diehl $'
#include "compilation_info.f90" #include "compilation_info.f90"
@ -121,26 +119,31 @@ subroutine AL_init()
! allocate (Fdot,source = F_lastInc) somethin like that should be possible ! allocate (Fdot,source = F_lastInc) somethin like that should be possible
allocate (F_lambda_lastInc(3,3, res(1), res(2),res(3)), source = 0.0_pReal) allocate (F_lambda_lastInc(3,3, res(1), res(2),res(3)), source = 0.0_pReal)
allocate (F_lambdaDot(3,3, res(1), res(2),res(3)), source = 0.0_pReal) allocate (F_lambdaDot(3,3, res(1), res(2),res(3)), source = 0.0_pReal)
allocate (P (3,3, res(1), res(2),res(3)), source = 0.0_pReal)
allocate (coordinates( res(1), res(2),res(3),3), source = 0.0_pReal) allocate (coordinates( res(1), res(2),res(3),3), source = 0.0_pReal)
allocate (temperature( res(1), res(2),res(3)), source = 0.0_pReal)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! PETSc Init ! PETSc Init
call SNESCreate(PETSC_COMM_WORLD,snes,ierr) call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
CHKERRQ(ierr)
call DMDACreate3d(PETSC_COMM_WORLD, & call DMDACreate3d(PETSC_COMM_WORLD, &
DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE, & DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE, &
DMDA_STENCIL_BOX,res(1),res(2),res(3),PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE, & DMDA_STENCIL_BOX,res(1),res(2),res(3),PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE, &
18,1,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr) 18,1,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr)
CHKERRQ(ierr)
call DMCreateGlobalVector(da,solution_vec,ierr) call DMCreateGlobalVector(da,solution_vec,ierr)
CHKERRQ(ierr)
call DMDASetLocalFunction(da,AL_formResidual,ierr) call DMDASetLocalFunction(da,AL_formResidual,ierr)
CHKERRQ(ierr)
call SNESSetDM(snes,da,ierr) call SNESSetDM(snes,da,ierr)
CHKERRQ(ierr)
call SNESSetConvergenceTest(snes,AL_converged,dummy,PETSC_NULL_FUNCTION,ierr) call SNESSetConvergenceTest(snes,AL_converged,dummy,PETSC_NULL_FUNCTION,ierr)
CHKERRQ(ierr)
call SNESSetFromOptions(snes,ierr) call SNESSetFromOptions(snes,ierr)
CHKERRQ(ierr)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! init fields ! init fields
call DMDAVecGetArrayF90(da,solution_vec,xx_psc,ierr) call DMDAVecGetArrayF90(da,solution_vec,xx_psc,ierr) ! places pointer xx_psc on PETSc data
F => xx_psc(0:8,:,:,:) F => xx_psc(0:8,:,:,:)
F_lambda => xx_psc(9:17,:,:,:) F_lambda => xx_psc(9:17,:,:,:)
if (restartInc == 1_pInt) then ! no deformation (no restart) if (restartInc == 1_pInt) then ! no deformation (no restart)
@ -177,12 +180,11 @@ subroutine AL_init()
call IO_read_jobBinaryFile(777,'F_aim_lastInc',trim(getSolverJobName()),size(F_aim_lastInc)) call IO_read_jobBinaryFile(777,'F_aim_lastInc',trim(getSolverJobName()),size(F_aim_lastInc))
read (777,rec=1) F_aim_lastInc read (777,rec=1) F_aim_lastInc
close (777) close (777)
coordinates = 0.0 ! change it later!!! coordinates = 0.0 ! change it later!!!
endif endif
!print*, shape(F)
call Utilities_constitutiveResponse(coordinates,F,F,temperature,0.0_pReal,P,C,P_av,.false.,math_I3) call Utilities_constitutiveResponse(coordinates,F,F,temperature,0.0_pReal,P,C,P_av,.false.,math_I3)
call DMDAVecRestoreArrayF90(da,solution_vec,xx_psc,ierr) call DMDAVecRestoreArrayF90(da,solution_vec,xx_psc,ierr)
CHKERRQ(ierr)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! reference stiffness ! reference stiffness
@ -202,6 +204,7 @@ subroutine AL_init()
end subroutine AL_init end subroutine AL_init
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief solution for the AL scheme with internal iterations !> @brief solution for the AL scheme with internal iterations
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
@ -226,7 +229,6 @@ subroutine AL_init()
Utilities_maskedCompliance, & Utilities_maskedCompliance, &
Utilities_updateGamma, & Utilities_updateGamma, &
cutBack cutBack
use FEsolving, only: & use FEsolving, only: &
restartWrite, & restartWrite, &
terminallyIll terminallyIll
@ -243,12 +245,13 @@ subroutine AL_init()
real(pReal), dimension(3,3), intent(in) :: rotation_BC real(pReal), dimension(3,3), intent(in) :: rotation_BC
real(pReal), dimension(3,3) ,save :: F_aimDot real(pReal), dimension(3,3) ,save :: F_aimDot
real(pReal), dimension(3,3) :: F_aim_lab real(pReal), dimension(3,3) :: F_aim_lab
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! loop variables, convergence etc. ! loop variables, convergence etc.
real(pReal), dimension(3,3) :: temp33_Real real(pReal), dimension(3,3) :: temp33_Real
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! ! PETSc Data
PetscScalar, dimension(:,:,:,:), pointer :: xx_psc, F, F_lambda PetscScalar, dimension(:,:,:,:), pointer :: xx_psc, F, F_lambda
PetscErrorCode :: ierr PetscErrorCode :: ierr
SNESConvergedReason ::reason SNESConvergedReason ::reason
@ -270,15 +273,12 @@ subroutine AL_init()
F => xx_psc(0:8,:,:,:) F => xx_psc(0:8,:,:,:)
F_lambda => xx_psc(9:17,:,:,:) F_lambda => xx_psc(9:17,:,:,:)
if ( cutBack) then if ( cutBack) then
F_aim = F_aim_lastInc F_aim = F_aim_lastInc
F_lambda= reshape(F_lambda_lastInc,[9,res(1),res(2),res(3)]) F_lambda= reshape(F_lambda_lastInc,[9,res(1),res(2),res(3)])
F = reshape(F_lastInc,[9,res(1),res(2),res(3)]) F = reshape(F_lastInc,[9,res(1),res(2),res(3)])
C = C_lastInc C = C_lastInc
else else
!--------------------------------------------------------------------------------------------------
C_lastInc = C C_lastInc = C
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! calculate rate for aim ! calculate rate for aim
@ -309,11 +309,11 @@ else
! update local deformation gradient and coordinates ! update local deformation gradient and coordinates
! deltaF_aim = math_rotate_backward33(deltaF_aim,rotation_BC) ! deltaF_aim = math_rotate_backward33(deltaF_aim,rotation_BC)
F = reshape(Utilities_forwardField(timeinc,F_aim,F_lastInc,Fdot),[9,res(1),res(2),res(3)]) F = reshape(Utilities_forwardField(timeinc,F_aim,F_lastInc,Fdot),[9,res(1),res(2),res(3)])
F_lambda = reshape(Utilities_forwardField(timeinc,F_aim,F_lambda_lastInc,F_lambdadot),[9,res(1),res(2),res(3)]) F_lambda = reshape(Utilities_forwardField(timeinc,F_aim,F_lambda_lastInc,F_lambdadot),[9,res(1),res(2),res(3)])
call DMDAVecRestoreArrayF90(da,solution_vec,xx_psc,ierr) call DMDAVecRestoreArrayF90(da,solution_vec,xx_psc,ierr)
CHKERRQ(ierr)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! update stiffness (and gamma operator) ! update stiffness (and gamma operator)
@ -327,18 +327,20 @@ else
params%timeinc = timeinc params%timeinc = timeinc
call SNESSolve(snes,PETSC_NULL_OBJECT,solution_vec,ierr) call SNESSolve(snes,PETSC_NULL_OBJECT,solution_vec,ierr)
CHKERRQ(ierr)
call SNESGetConvergedReason(snes,reason,ierr) call SNESGetConvergedReason(snes,reason,ierr)
CHKERRQ(ierr)
AL_solution%termIll = terminallyIll AL_solution%termIll = terminallyIll
terminallyIll = .false. terminallyIll = .false.
if (reason > 0 ) then if (reason > 0 ) then
AL_solution%converged = .true. AL_solution%converged = .true.
AL_solution%iterationsNeeded = reportIter AL_solution%iterationsNeeded = reportIter
endif endif
end function AL_solution end function AL_solution
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief forms the AL residual vector !> @brief forms the AL residual vector
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
@ -363,7 +365,6 @@ else
use IO, only: IO_intOut use IO, only: IO_intOut
implicit none implicit none
integer(pInt) :: i,j,k integer(pInt) :: i,j,k
integer(pInt), save :: callNo = 3_pInt integer(pInt), save :: callNo = 3_pInt
real(pReal), dimension(3,3) :: temp33_Real real(pReal), dimension(3,3) :: temp33_Real
@ -384,7 +385,9 @@ else
residual_F_lambda => f_scal(:,:,2,:,:,:) residual_F_lambda => f_scal(:,:,2,:,:,:)
call SNESGetNumberFunctionEvals(snes,nfuncs,ierr) call SNESGetNumberFunctionEvals(snes,nfuncs,ierr)
CHKERRQ(ierr)
call SNESGetIterationNumber(snes,iter,ierr) call SNESGetIterationNumber(snes,iter,ierr)
CHKERRQ(ierr)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! report begin of new iteration ! report begin of new iteration
@ -400,9 +403,9 @@ else
reportIter = reportIter + 1_pInt reportIter = reportIter + 1_pInt
endif endif
callNo = callNo +1_pInt callNo = callNo +1_pInt
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! evaluate constitutive response ! evaluate constitutive response
call Utilities_constitutiveResponse(coordinates,F_lastInc,F,temperature,params%timeinc, & call Utilities_constitutiveResponse(coordinates,F_lastInc,F,temperature,params%timeinc, &
residual_F,C,P_av,ForwardData,params%rotation_BC) residual_F,C,P_av,ForwardData,params%rotation_BC)
ForwardData = .False. ForwardData = .False.
@ -412,9 +415,8 @@ else
F_aim = F_aim - math_mul3333xx33(S, ((P_av - params%P_BC))) ! S = 0.0 for no bc F_aim = F_aim - math_mul3333xx33(S, ((P_av - params%P_BC))) ! S = 0.0 for no bc
err_stress = maxval(abs(mask_stress * (P_av - params%P_BC))) ! mask = 0.0 for no bc err_stress = maxval(abs(mask_stress * (P_av - params%P_BC))) ! mask = 0.0 for no bc
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! doing Fourier transform !
field_real = 0.0_pReal field_real = 0.0_pReal
do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1) do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1)
temp33_Real = math_mul3333xx33(S_scale,residual_F(1:3,1:3,i,j,k)) + math_I3 temp33_Real = math_mul3333xx33(S_scale,residual_F(1:3,1:3,i,j,k)) + math_I3
@ -423,7 +425,7 @@ else
enddo; enddo; enddo enddo; enddo; enddo
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! doing Fourier transform ! doing convolution in Fourier space
call Utilities_FFTforward() call Utilities_FFTforward()
call Utilities_fourierConvolution(math_rotate_backward33(F_aim,params%rotation_BC)) call Utilities_fourierConvolution(math_rotate_backward33(F_aim,params%rotation_BC))
call Utilities_FFTbackward() call Utilities_FFTbackward()
@ -441,6 +443,7 @@ else
end subroutine AL_formResidual end subroutine AL_formResidual
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief convergence check !> @brief convergence check
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
@ -477,9 +480,12 @@ else
reason = 0 reason = 0
endif endif
write(6,'(a,es14.7)') 'error stress BC = ', err_stress/min(maxval(abs(P_av))*err_stress_tolrel,err_stress_tolabs) write(6,'(a,es14.7)') 'error stress BC = ', &
write(6,'(a,es14.7)') 'error F = ', err_f/sqrt(sum((F_aim-math_I3)*(F_aim-math_I3)))/err_f_tol err_stress/min(maxval(abs(P_av))*err_stress_tolrel,err_stress_tolabs)
write(6,'(a,es14.7)') 'error P = ', err_p/sqrt(sum((F_aim-math_I3)*(F_aim-math_I3)))/err_p_tol write(6,'(a,es14.7)') 'error F = ',&
err_f/sqrt(sum((F_aim-math_I3)*(F_aim-math_I3)))/err_f_tol
write(6,'(a,es14.7)') 'error P = ', &
err_p/sqrt(sum((F_aim-math_I3)*(F_aim-math_I3)))/err_p_tol
write(6,'(/,a)') '==========================================================================' write(6,'(/,a)') '=========================================================================='
end subroutine AL_converged end subroutine AL_converged
@ -494,9 +500,13 @@ else
PetscErrorCode :: ierr PetscErrorCode :: ierr
call VecDestroy(solution_vec,ierr) call VecDestroy(solution_vec,ierr)
CHKERRQ(ierr)
call SNESDestroy(snes,ierr) call SNESDestroy(snes,ierr)
CHKERRQ(ierr)
call DMDestroy(da,ierr) call DMDestroy(da,ierr)
CHKERRQ(ierr)
call PetscFinalize(ierr) call PetscFinalize(ierr)
CHKERRQ(ierr)
call Utilities_destroy() call Utilities_destroy()
end subroutine AL_destroy end subroutine AL_destroy

View File

@ -12,10 +12,8 @@ module DAMASK_spectral_SolverBasic
use prec, only: & use prec, only: &
pInt, & pInt, &
pReal pReal
use math, only: & use math, only: &
math_I3 math_I3
use DAMASK_spectral_Utilities, only: & use DAMASK_spectral_Utilities, only: &
tSolutionState tSolutionState
@ -24,19 +22,23 @@ module DAMASK_spectral_SolverBasic
DAMASK_spectral_SolverBasic_label = 'basic' DAMASK_spectral_SolverBasic_label = 'basic'
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! pointwise data ! pointwise global data
real(pReal), private, dimension(:,:,:,:,:), allocatable :: F, F_lastInc, Fdot, P real(pReal), private, dimension(:,:,:,:,:), allocatable :: &
F, & !< deformation gradient field
F_lastInc, & !< deformation gradient field last increment
Fdot !< assumed rate for F n to F n+1
real(pReal), private, dimension(:,:,:,:), allocatable :: coordinates real(pReal), private, dimension(:,:,:,:), allocatable :: coordinates
real(pReal), private, dimension(:,:,:), allocatable :: temperature real(pReal), private :: temperature !< not pointwise
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! stress, stiffness and compliance average etc. ! stress, stiffness and compliance average etc.
real(pReal), private, dimension(3,3) :: & real(pReal), private, dimension(3,3) :: &
F_aim = math_I3, & F_aim = math_I3, & !< deformation gradient aim
F_aim_lastInc = math_I3, & F_aim_lastInc = math_I3, & !< deformation gradient aim last increment
F_aimDot = 0.0_pReal F_aimDot = 0.0_pReal !< assumed rate
real(pReal), private,dimension(3,3,3,3) :: & real(pReal), private,dimension(3,3,3,3) :: &
C = 0.0_pReal, C_lastInc = 0.0_pReal C = 0.0_pReal, & !< average stiffness
C_lastInc = 0.0_pReal !< average stiffness last increment
contains contains
@ -50,27 +52,27 @@ subroutine basic_init()
IO_read_JobBinaryFile, & IO_read_JobBinaryFile, &
IO_write_JobBinaryFile, & IO_write_JobBinaryFile, &
IO_intOut IO_intOut
use FEsolving, only: & use FEsolving, only: &
restartInc restartInc
use DAMASK_interface, only: & use DAMASK_interface, only: &
getSolverJobName getSolverJobName
use DAMASK_spectral_Utilities, only: & use DAMASK_spectral_Utilities, only: &
Utilities_init, & Utilities_init, &
Utilities_constitutiveResponse, & Utilities_constitutiveResponse, &
Utilities_updateGamma, & Utilities_updateGamma, &
debugRestart debugRestart
use mesh, only: & use mesh, only: &
res, & res, &
geomdim geomdim
implicit none implicit none
integer(pInt) :: i,j,k real(pReal), dimension(3,3,res(1),res(2),res(3)) :: P
real(pReal), dimension(3,3) :: temp33_Real integer(pInt) :: &
real(pReal), dimension(3,3,3,3) :: temp3333_Real i, j, k
real(pReal), dimension(3,3) :: &
temp33_Real
real(pReal), dimension(3,3,3,3) :: &
temp3333_Real
call Utilities_Init() call Utilities_Init()
@ -79,15 +81,15 @@ subroutine basic_init()
#include "compilation_info.f90" #include "compilation_info.f90"
write(6,'(a)') '' write(6,'(a)') ''
!--------------------------------------------------------------------------------------------------
! allocate global fields
allocate (F ( 3,3,res(1), res(2),res(3)), source = 0.0_pReal) allocate (F ( 3,3,res(1), res(2),res(3)), source = 0.0_pReal)
allocate (F_lastInc ( 3,3,res(1), res(2),res(3)), source = 0.0_pReal) allocate (F_lastInc ( 3,3,res(1), res(2),res(3)), source = 0.0_pReal)
allocate (Fdot ( 3,3,res(1), res(2),res(3)), source = 0.0_pReal) allocate (Fdot ( 3,3,res(1), res(2),res(3)), source = 0.0_pReal)
allocate (P ( 3,3,res(1), res(2),res(3)), source = 0.0_pReal)
allocate (coordinates( res(1), res(2),res(3),3),source = 0.0_pReal) allocate (coordinates( res(1), res(2),res(3),3),source = 0.0_pReal)
allocate (temperature( res(1), res(2),res(3)), source = 0.0_pReal)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! init fields ! init fields and average quantities
if (restartInc == 1_pInt) then ! no deformation (no restart) if (restartInc == 1_pInt) then ! no deformation (no restart)
F = spread(spread(spread(math_I3,3,res(1)),4,res(2)),5,res(3)) ! initialize to identity F = spread(spread(spread(math_I3,3,res(1)),4,res(2)),5,res(3)) ! initialize to identity
F_lastInc = F F_lastInc = F
@ -127,13 +129,9 @@ subroutine basic_init()
close (777) close (777)
coordinates = 0.0 ! change it later!!! coordinates = 0.0 ! change it later!!!
endif endif
call Utilities_constitutiveResponse(coordinates,F,F,temperature,0.0_pReal,& call Utilities_constitutiveResponse(coordinates,F,F,temperature,0.0_pReal,P,C,temp33_Real,.false.,math_I3) ! constitutive response with no deformation in no time to get reference stiffness
P,C,temp33_Real,.false.,math_I3)
!no rotation bc call deformed_fft(res,geomdim,math_rotate_backward33(F_aim,rotation_BC),1.0_pReal,F_lastInc,coordinates) !no rotation bc call deformed_fft(res,geomdim,math_rotate_backward33(F_aim,rotation_BC),1.0_pReal,F_lastInc,coordinates)
!--------------------------------------------------------------------------------------------------
! reference stiffness
if (restartInc == 1_pInt) then ! use initial stiffness as reference stiffness if (restartInc == 1_pInt) then ! use initial stiffness as reference stiffness
temp3333_Real = C temp3333_Real = C
endif endif
@ -148,7 +146,6 @@ end subroutine basic_init
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
type(tSolutionState) function & type(tSolutionState) function &
basic_solution(incInfo,guess,timeinc,timeinc_old,P_BC,F_BC,temperature_bc,rotation_BC) basic_solution(incInfo,guess,timeinc,timeinc_old,P_BC,F_BC,temperature_bc,rotation_BC)
use numerics, only: & use numerics, only: &
itmax, & itmax, &
itmin, & itmin, &
@ -166,7 +163,6 @@ type(tSolutionState) function &
use IO, only: & use IO, only: &
IO_write_JobBinaryFile, & IO_write_JobBinaryFile, &
IO_intOut IO_intOut
use DAMASK_spectral_Utilities, only: & use DAMASK_spectral_Utilities, only: &
tBoundaryCondition, & tBoundaryCondition, &
field_real, & field_real, &
@ -179,7 +175,6 @@ type(tSolutionState) function &
Utilities_updateGamma, & Utilities_updateGamma, &
Utilities_constitutiveResponse, & Utilities_constitutiveResponse, &
Utilities_calculateRate Utilities_calculateRate
use FEsolving, only: & use FEsolving, only: &
restartWrite, & restartWrite, &
restartRead, & restartRead, &
@ -190,26 +185,40 @@ type(tSolutionState) function &
implicit none implicit none
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! input data for solution ! input data for solution
real(pReal), intent(in) :: timeinc, timeinc_old, temperature_bc real(pReal), intent(in) :: &
logical, intent(in) :: guess timeinc, & !< increment in time for current solution
type(tBoundaryCondition), intent(in) :: P_BC,F_BC timeinc_old, & !< increment in time of last increment
character(len=*), intent(in) :: incInfo temperature_bc !< temperature (I wonder, why it is intent(in)
real(pReal), dimension(3,3), intent(in) :: rotation_BC logical, intent(in) :: &
guess !< if .false., assume homogeneous addon
type(tBoundaryCondition), intent(in) :: &
P_BC,& !< stress boundary conditions
F_BC !< deformation boundary conditions
character(len=*), intent(in) :: &
incInfo !< string with information of current increment for output to screen
real(pReal), dimension(3,3), intent(in) :: &
rotation_BC !< rotation load to lab
real(pReal), dimension(3,3,3,3) :: S real(pReal), dimension(3,3,3,3) :: &
real(pReal), dimension(3,3) :: F_aim_lab, & S !< current average compliance
F_aim_lab_lastIter, & real(pReal), dimension(3,3) :: &
F_aim_lab, &
F_aim_lab_lastIter, & !< aim of last iteration
P_av P_av
real(pReal), dimension(3,3,res(1),res(2),res(3)) :: P
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! loop variables, convergence etc. ! loop variables, convergence etc.
real(pReal) :: err_div, err_stress real(pReal) :: err_div, err_stress
integer(pInt) :: iter, row, column, i, j, k integer(pInt) :: iter, row, column, i, j, k
logical :: ForwardData logical :: ForwardData
real(pReal) :: defgradDet, defgradDetMax, defgradDetMin real(pReal) :: &
defgradDet, &
defgradDetMax, &
defgradDetMin
real(pReal), dimension(3,3) :: temp33_Real real(pReal), dimension(3,3) :: temp33_Real
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! restart information for spectral solver ! write restart information for spectral solver
if (restartWrite) then if (restartWrite) then
write(6,'(a)') 'writing converged results for restart' write(6,'(a)') 'writing converged results for restart'
call IO_write_jobBinaryFile(777,'convergedSpectralDefgrad',size(F)) ! writing deformation gradient field to file call IO_write_jobBinaryFile(777,'convergedSpectralDefgrad',size(F)) ! writing deformation gradient field to file
@ -235,14 +244,14 @@ type(tSolutionState) function &
close(777) close(777)
endif endif
!--------------------------------------------------------------------------------------------------
! forward data
if (cutBack) then if (cutBack) then
F_aim = F_aim_lastInc F_aim = F_aim_lastInc
F = F_lastInc F = F_lastInc
C = C_lastInc C = C_lastInc
else else
C_lastInc = C C_lastInc = C
!--------------------------------------------------------------------------------------------------
! calculate rate for aim
if (F_BC%myType=='l') then ! calculate f_aimDot from given L and current F if (F_BC%myType=='l') then ! calculate f_aimDot from given L and current F
f_aimDot = F_BC%maskFloat * math_mul33x33(F_BC%values, F_aim) f_aimDot = F_BC%maskFloat * math_mul33x33(F_BC%values, F_aim)
elseif(F_BC%myType=='fdot') then ! f_aimDot is prescribed elseif(F_BC%myType=='fdot') then ! f_aimDot is prescribed
@ -251,9 +260,7 @@ type(tSolutionState) function &
if (guess) f_aimDot = f_aimDot + P_BC%maskFloat * (F_aim - F_aim_lastInc)/timeinc_old if (guess) f_aimDot = f_aimDot + P_BC%maskFloat * (F_aim - F_aim_lastInc)/timeinc_old
F_aim_lastInc = F_aim F_aim_lastInc = F_aim
!-------------------------------------------------------------------------------------------------- call deformed_fft(res,geomdim,math_rotate_backward33(F_aim_lastInc,rotation_BC), & ! update coordinates and rate and forward last inc
! update coordinates and rate and forward last inc
call deformed_fft(res,geomdim,math_rotate_backward33(F_aim_lastInc,rotation_BC), &
1.0_pReal,F_lastInc,coordinates) 1.0_pReal,F_lastInc,coordinates)
Fdot = Utilities_calculateRate(math_rotate_backward33(f_aimDot,rotation_BC), & Fdot = Utilities_calculateRate(math_rotate_backward33(f_aimDot,rotation_BC), &
timeinc,timeinc_old,guess,F_lastInc,F) timeinc,timeinc_old,guess,F_lastInc,F)
@ -265,13 +272,13 @@ type(tSolutionState) function &
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! update stiffness (and gamma operator) ! update stiffness (and gamma operator)
S = Utilities_maskedCompliance(rotation_BC,P_BC%maskLogical,C) S = Utilities_maskedCompliance(rotation_BC,P_BC%maskLogical,C)
if (update_gamma) call Utilities_updateGamma(C,restartWrite) if (update_gamma) call Utilities_updateGamma(C,restartWrite)
!--------------------------------------------------------------------------------------------------
! iteration till converged
if (.not. restartRead) ForwardData = .True. if (.not. restartRead) ForwardData = .True.
iter = 0_pInt iter = 0_pInt
convergenceLoop: do while(iter < itmax) convergenceLoop: do while(iter < itmax)
iter = iter + 1_pInt iter = iter + 1_pInt
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! report begin of new iteration ! report begin of new iteration
@ -279,16 +286,16 @@ type(tSolutionState) function &
' @ Iter. ', itmin, '≤',iter, '≤', itmax ' @ Iter. ', itmin, '≤',iter, '≤', itmax
write(6,'(a,/,3(3(f12.7,1x)/))',advance='no') 'deformation gradient aim =', & write(6,'(a,/,3(3(f12.7,1x)/))',advance='no') 'deformation gradient aim =', &
math_transpose33(F_aim) math_transpose33(F_aim)
F_aim_lab_lastIter = math_rotate_backward33(F_aim,rotation_BC)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! evaluate constitutive response ! evaluate constitutive response
F_aim_lab_lastIter = math_rotate_backward33(F_aim,rotation_BC)
basic_solution%termIll = .false. basic_solution%termIll = .false.
call Utilities_constitutiveResponse(coordinates,F_lastInc,F,temperature,timeinc,& call Utilities_constitutiveResponse(coordinates,F_lastInc,F,temperature,timeinc,&
P,C,P_av,ForwardData,rotation_BC) P,C,P_av,ForwardData,rotation_BC)
basic_solution%termIll = terminallyIll basic_solution%termIll = terminallyIll
terminallyIll = .False. terminallyIll = .false.
ForwardData = .False. ForwardData = .false.
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! stress BC handling ! stress BC handling
@ -326,7 +333,6 @@ logical function basic_Converged(err_div,pAvgDiv,err_stress,pAvgStress)
err_div_tol, & err_div_tol, &
err_stress_tolrel, & err_stress_tolrel, &
err_stress_tolabs err_stress_tolabs
use math, only: & use math, only: &
math_mul33x33, & math_mul33x33, &
math_eigenvalues33, & math_eigenvalues33, &
@ -345,8 +351,6 @@ logical function basic_Converged(err_div,pAvgDiv,err_stress,pAvgStress)
err_stress_tol, & err_stress_tol, &
pAvgDivL2 pAvgDivL2
pAvgDivL2 = sqrt(maxval(math_eigenvalues33(math_mul33x33(pAvgDiv,math_transpose33(pAvgDiv))))) ! L_2 norm of average stress (http://mathworld.wolfram.com/SpectralNorm.html) pAvgDivL2 = sqrt(maxval(math_eigenvalues33(math_mul33x33(pAvgDiv,math_transpose33(pAvgDiv))))) ! L_2 norm of average stress (http://mathworld.wolfram.com/SpectralNorm.html)
err_stress_tol = min(maxval(abs(pAvgStress))*err_stress_tolrel,err_stress_tolabs) err_stress_tol = min(maxval(abs(pAvgStress))*err_stress_tolrel,err_stress_tolabs)
@ -372,44 +376,4 @@ subroutine basic_destroy()
end subroutine basic_destroy end subroutine basic_destroy
end module DAMASK_spectral_SolverBasic end module DAMASK_spectral_SolverBasic
!--------------------------------------------------------------------------------------------------
! calculate some additional output
! if(debugGeneral) then
! maxCorrectionSkew = 0.0_pReal
! maxCorrectionSym = 0.0_pReal
! temp33_Real = 0.0_pReal
! do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1)
! maxCorrectionSym = max(maxCorrectionSym,&
! maxval(math_symmetric33(field_real(i,j,k,1:3,1:3))))
! maxCorrectionSkew = max(maxCorrectionSkew,&
! maxval(math_skew33(field_real(i,j,k,1:3,1:3))))
! temp33_Real = temp33_Real + field_real(i,j,k,1:3,1:3)
! enddo; enddo; enddo
! write(6,'(a,1x,es11.4)') 'max symmetric correction of deformation =',&
! maxCorrectionSym*wgt
! write(6,'(a,1x,es11.4)') 'max skew correction of deformation =',&
! maxCorrectionSkew*wgt
! write(6,'(a,1x,es11.4)') 'max sym/skew of avg correction = ',&
! maxval(math_symmetric33(temp33_real))/&
! maxval(math_skew33(temp33_real))
! endif
!--------------------------------------------------------------------------------------------------
! calculate bounds of det(F) and report
! if(debugGeneral) then
! defgradDetMax = -huge(1.0_pReal)
! defgradDetMin = +huge(1.0_pReal)
! do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1)
! defgradDet = math_det33(F(i,j,k,1:3,1:3))
! defgradDetMax = max(defgradDetMax,defgradDet)
! defgradDetMin = min(defgradDetMin,defgradDet)
! enddo; enddo; enddo
! write(6,'(a,1x,es11.4)') 'max determinant of deformation =', defgradDetMax
! write(6,'(a,1x,es11.4)') 'min determinant of deformation =', defgradDetMin
! endif

View File

@ -45,7 +45,7 @@ module DAMASK_spectral_SolverBasicPETSc
! common pointwise data ! common pointwise data
real(pReal), private, dimension(:,:,:,:,:), allocatable :: F_lastInc, Fdot real(pReal), private, dimension(:,:,:,:,:), allocatable :: F_lastInc, Fdot
real(pReal), private, dimension(:,:,:,:), allocatable :: coordinates real(pReal), private, dimension(:,:,:,:), allocatable :: coordinates
real(pReal), private, dimension(:,:,:), allocatable :: temperature real(pReal) :: temperature
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! stress, stiffness and compliance average etc. ! stress, stiffness and compliance average etc.
@ -67,6 +67,7 @@ module DAMASK_spectral_SolverBasicPETSc
basicPETSc_solution ,& basicPETSc_solution ,&
basicPETSc_destroy basicPETSc_destroy
contains contains
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief allocates all neccessary fields and fills them with data, potentially from restart info !> @brief allocates all neccessary fields and fills them with data, potentially from restart info
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
@ -111,17 +112,19 @@ subroutine basicPETSc_init()
PetscObject :: dummy PetscObject :: dummy
call Utilities_init() call Utilities_init()
write(6,'(/,a)') ' <<<+- DAMASK_spectral_solverBasicPETSc init -+>>>' write(6,'(/,a)') ' <<<+- DAMASK_spectral_solverBasicPETSc init -+>>>'
write(6,'(a)') ' $Id: DAMASK_spectral_SolverBasicPETSC.f90 1654 2012-08-03 09:25:48Z MPIE\m.diehl $' write(6,'(a)') ' $Id: DAMASK_spectral_SolverBasicPETSC.f90 1654 2012-08-03 09:25:48Z MPIE\m.diehl $'
#include "compilation_info.f90" #include "compilation_info.f90"
write(6,'(a)') '' write(6,'(a)') ''
!--------------------------------------------------------------------------------------------------
! allocate global fields
allocate (F_lastInc (3,3, res(1), res(2),res(3)), source = 0.0_pReal) allocate (F_lastInc (3,3, res(1), res(2),res(3)), source = 0.0_pReal)
allocate (Fdot (3,3, res(1), res(2),res(3)), source = 0.0_pReal) allocate (Fdot (3,3, res(1), res(2),res(3)), source = 0.0_pReal)
allocate (coordinates( res(1), res(2),res(3),3), source = 0.0_pReal) allocate (coordinates( res(1), res(2),res(3),3), source = 0.0_pReal)
allocate (temperature( res(1), res(2),res(3)), source = 0.0_pReal)
!--------------------------------------------------------------------------------------------------
! initialize solver specific parts of PETSc
call SNESCreate(PETSC_COMM_WORLD,snes,ierr) call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
CHKERRQ(ierr) CHKERRQ(ierr)
call DMDACreate3d(PETSC_COMM_WORLD, & call DMDACreate3d(PETSC_COMM_WORLD, &
@ -142,7 +145,7 @@ subroutine basicPETSc_init()
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! init fields ! init fields
call DMDAVecGetArrayF90(da,solution_vec,F,ierr) call DMDAVecGetArrayF90(da,solution_vec,F,ierr) ! get the data out of PETSc to work with
CHKERRQ(ierr) CHKERRQ(ierr)
if (restartInc == 1_pInt) then ! no deformation (no restart) if (restartInc == 1_pInt) then ! no deformation (no restart)
F_lastInc = spread(spread(spread(math_I3,3,res(1)),4,res(2)),5,res(3)) ! initialize to identity F_lastInc = spread(spread(spread(math_I3,3,res(1)),4,res(2)),5,res(3)) ! initialize to identity
@ -177,12 +180,11 @@ subroutine basicPETSc_init()
reshape(F(0:8,0:res(1)-1_pInt,0:res(2)-1_pInt,0:res(3)-1_pInt),[3,3,res(1),res(2),res(3)]),& reshape(F(0:8,0:res(1)-1_pInt,0:res(2)-1_pInt,0:res(3)-1_pInt),[3,3,res(1),res(2),res(3)]),&
reshape(F(0:8,0:res(1)-1_pInt,0:res(2)-1_pInt,0:res(3)-1_pInt),[3,3,res(1),res(2),res(3)]),& reshape(F(0:8,0:res(1)-1_pInt,0:res(2)-1_pInt,0:res(3)-1_pInt),[3,3,res(1),res(2),res(3)]),&
temperature,0.0_pReal,P,C,P_av,.false.,math_I3) temperature,0.0_pReal,P,C,P_av,.false.,math_I3)
call DMDAVecRestoreArrayF90(da,solution_vec,F,ierr) ! write data back into PETSc
call DMDAVecRestoreArrayF90(da,solution_vec,F,ierr) CHKERRQ(ierr)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! reference stiffness ! reference stiffness and Gamma operator
if (restartInc == 1_pInt) then if (restartInc == 1_pInt) then
call IO_write_jobBinaryFile(777,'C_ref',size(C)) call IO_write_jobBinaryFile(777,'C_ref',size(C))
write (777,rec=1) C write (777,rec=1) C
@ -192,7 +194,6 @@ subroutine basicPETSc_init()
read (777,rec=1) C read (777,rec=1) C
close (777) close (777)
endif endif
call Utilities_updateGamma(C,.True.) call Utilities_updateGamma(C,.True.)
end subroutine basicPETSc_init end subroutine basicPETSc_init
@ -265,7 +266,6 @@ if ( cutBack) then
F = reshape(F_lastInc,[9,res(1),res(2),res(3)]) F = reshape(F_lastInc,[9,res(1),res(2),res(3)])
C = C_lastInc C = C_lastInc
else else
C_lastInc = C C_lastInc = C
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
@ -291,6 +291,7 @@ else
F = reshape(Utilities_forwardField(timeinc,F_aim,F_lastInc,Fdot),[9,res(1),res(2),res(3)]) F = reshape(Utilities_forwardField(timeinc,F_aim,F_lastInc,Fdot),[9,res(1),res(2),res(3)])
call DMDAVecRestoreArrayF90(da,solution_vec,F,ierr) call DMDAVecRestoreArrayF90(da,solution_vec,F,ierr)
CHKERRQ(ierr)
call deformed_fft(res,geomdim,math_rotate_backward33(F_aim,rotation_BC),1.0_pReal,F_lastInc,coordinates) call deformed_fft(res,geomdim,math_rotate_backward33(F_aim,rotation_BC),1.0_pReal,F_lastInc,coordinates)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
@ -305,7 +306,9 @@ else
params%timeinc = timeinc params%timeinc = timeinc
call SNESSolve(snes,PETSC_NULL_OBJECT,solution_vec,ierr) call SNESSolve(snes,PETSC_NULL_OBJECT,solution_vec,ierr)
CHKERRQ(ierr)
call SNESGetConvergedReason(snes,reason,ierr) call SNESGetConvergedReason(snes,reason,ierr)
CHKERRQ(ierr)
basicPETSc_solution%termIll = terminallyIll basicPETSc_solution%termIll = terminallyIll
terminallyIll = .false. terminallyIll = .false.
BasicPETSC_solution%converged =.false. BasicPETSC_solution%converged =.false.
@ -338,8 +341,8 @@ else
Utilities_constitutiveResponse, & Utilities_constitutiveResponse, &
Utilities_divergenceRMS Utilities_divergenceRMS
use IO, only : IO_intOut use IO, only : IO_intOut
implicit none
implicit none
real(pReal), dimension(3,3) :: F_aim_lab_lastIter, F_aim_lab real(pReal), dimension(3,3) :: F_aim_lab_lastIter, F_aim_lab
DMDALocalInfo, dimension(*) :: myIn DMDALocalInfo, dimension(*) :: myIn
@ -350,8 +353,9 @@ else
PetscErrorCode :: ierr PetscErrorCode :: ierr
call SNESGetNumberFunctionEvals(snes,nfuncs,ierr) call SNESGetNumberFunctionEvals(snes,nfuncs,ierr)
CHKERRQ(ierr)
call SNESGetIterationNumber(snes,iter,ierr) call SNESGetIterationNumber(snes,iter,ierr)
CHKERRQ(ierr)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! report begin of new iteration ! report begin of new iteration
write(6,'(/,a,3(a,'//IO_intOut(itmax)//'))') trim(incInfo), & write(6,'(/,a,3(a,'//IO_intOut(itmax)//'))') trim(incInfo), &
@ -388,6 +392,7 @@ else
write(6,'(/,a)') '==========================================================================' write(6,'(/,a)') '=========================================================================='
end subroutine BasicPETSc_formResidual end subroutine BasicPETSc_formResidual
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief convergence check !> @brief convergence check
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
@ -433,9 +438,9 @@ else
' (',err_div/pAvgDivL2,' N/m³)' ' (',err_div/pAvgDivL2,' N/m³)'
write(6,'(a,f6.2,a,es11.4,a)') 'error stress = ', err_stress/min(maxval(abs(P_av))*err_stress_tolrel,err_stress_tolabs), & write(6,'(a,f6.2,a,es11.4,a)') 'error stress = ', err_stress/min(maxval(abs(P_av))*err_stress_tolrel,err_stress_tolabs), &
' (',err_stress,' Pa)' ' (',err_stress,' Pa)'
end subroutine BasicPETSc_converged end subroutine BasicPETSc_converged
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief destroy routine !> @brief destroy routine
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
@ -443,15 +448,20 @@ else
use DAMASK_spectral_Utilities, only: & use DAMASK_spectral_Utilities, only: &
Utilities_destroy Utilities_destroy
implicit none implicit none
PetscErrorCode :: ierr PetscErrorCode :: ierr
call VecDestroy(solution_vec,ierr) call VecDestroy(solution_vec,ierr)
CHKERRQ(ierr)
call SNESDestroy(snes,ierr) call SNESDestroy(snes,ierr)
CHKERRQ(ierr)
call DMDestroy(da,ierr) call DMDestroy(da,ierr)
CHKERRQ(ierr)
call PetscFinalize(ierr) call PetscFinalize(ierr)
CHKERRQ(ierr)
call Utilities_destroy() call Utilities_destroy()
CHKERRQ(ierr)
end subroutine BasicPETSc_destroy end subroutine BasicPETSc_destroy
end module DAMASK_spectral_SolverBasicPETSc end module DAMASK_spectral_SolverBasicPETSc

View File

@ -17,26 +17,29 @@ module DAMASK_spectral_utilities
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! variables storing information for spectral method and FFTW ! variables storing information for spectral method and FFTW
real(pReal), public, dimension(:,:,:,:,:), pointer :: field_real !< real representation (some stress of deformation) of field_fourier real(pReal), public, dimension(:,:,:,:,:), pointer :: field_real !< real representation (some stress or deformation) of field_fourier
type(C_PTR), private :: plan_forward, plan_backward !< plans for FFTW complex(pReal),private, dimension(:,:,:,:,:), pointer :: field_fourier !< field on which the Fourier transform operates
real(pReal), private, dimension(:,:,:,:,:,:,:), allocatable :: gamma_hat !< gamma operator (field) for spectral method real(pReal), private, dimension(:,:,:,:,:,:,:), allocatable :: gamma_hat !< gamma operator (field) for spectral method
real(pReal), private, dimension(:,:,:,:), allocatable :: xi !< wave vector field for divergence and for gamma operator real(pReal), private, dimension(:,:,:,:), allocatable :: xi !< wave vector field for divergence and for gamma operator
complex(pReal),private, dimension(:,:,:,:,:), pointer :: field_fourier !< field on which the Fourier transform operates
real(pReal), private, dimension(3,3,3,3) :: C_ref !< reference stiffness real(pReal), private, dimension(3,3,3,3) :: C_ref !< reference stiffness
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! debug fftw ! debug fftw
type(C_PTR), private :: plan_scalarField_forth, plan_scalarField_back !< plans for FFTW in case of debugging the Fourier transform complex(pReal),private, dimension(:,:,:), pointer :: scalarField_real, & !< scalar field real representation for debug of FFTW
complex(pReal),private, dimension(:,:,:), pointer :: scalarField_real !< scalar field real representation for debug of FFTW scalarField_fourier !< scalar field complex representation for debug of FFTW
complex(pReal),private, dimension(:,:,:), pointer :: scalarField_fourier !< scalar field complex representation for debug of FFTW
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! debug divergence ! debug divergence
type(C_PTR), private :: plan_divergence !< plan for FFTW in case of debugging divergence calculation
real(pReal), private, dimension(:,:,:,:), pointer :: divergence_real !< scalar field real representation for debugging divergence calculation real(pReal), private, dimension(:,:,:,:), pointer :: divergence_real !< scalar field real representation for debugging divergence calculation
complex(pReal),private, dimension(:,:,:,:), pointer :: divergence_fourier !< scalar field real representation for debugging divergence calculation complex(pReal),private, dimension(:,:,:,:), pointer :: divergence_fourier !< scalar field real representation for debugging divergence calculation
real(pReal), private, dimension(:,:,:,:), allocatable :: divergence_post !< data of divergence calculation using function from core modules (serves as a reference) real(pReal), private, dimension(:,:,:,:), allocatable :: divergence_post !< data of divergence calculation using function from core modules (serves as a reference)
!--------------------------------------------------------------------------------------------------
! plans for FFTW
type(C_PTR), private :: plan_scalarField_forth, plan_scalarField_back !< plans for FFTW in case of debugging the Fourier transform
type(C_PTR), private :: plan_forward, plan_backward !< plans for FFTW
type(C_PTR), private :: plan_divergence !< plan for FFTW in case of debugging divergence calculation
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! variables controlling debugging ! variables controlling debugging
logical, public :: & logical, public :: &
@ -80,8 +83,8 @@ module DAMASK_spectral_utilities
contains contains
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief allocates all neccessary fields, sets debug flags, create plans for fftw !> @brief allocates all neccessary fields, sets debug flags, create plans for FFTW
!> @details Sets the debug levels for general, divergence, restart and fftw from the biwise coding !> @details Sets the debug levels for general, divergence, restart and FFTW from the biwise coding
!> provided by the debug module to logicals. !> provided by the debug module to logicals.
!> Allocates all fields used by FFTW and create the corresponding plans depending on the debug !> Allocates all fields used by FFTW and create the corresponding plans depending on the debug
!> level chosen. !> level chosen.
@ -132,8 +135,8 @@ subroutine utilities_init()
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! allocation ! allocation
allocate (xi (3,res1_red,res(2),res(3)), source = 0.0_pReal) ! start out isothermally allocate (xi(3,res1_red,res(2),res(3)),source = 0.0_pReal) ! frequencies, only half the size for first dimension
tensorField = fftw_alloc_complex(int(res1_red*res(2)*res(3)*9_pInt,C_SIZE_T)) ! allocate continous data using a C function, C_SIZE_T is of type integer(8) tensorField = fftw_alloc_complex(int(res1_red*res(2)*res(3)*9_pInt,C_SIZE_T)) ! allocate aligned data using a C function, C_SIZE_T is of type integer(8)
call c_f_pointer(tensorField, field_real, [ res(1)+2_pInt,res(2),res(3),3,3]) ! place a pointer for a real representation on tensorField call c_f_pointer(tensorField, field_real, [ res(1)+2_pInt,res(2),res(3),3,3]) ! place a pointer for a real representation on tensorField
call c_f_pointer(tensorField, field_fourier, [ res1_red, res(2),res(3),3,3]) ! place a pointer for a complex representation on tensorField call c_f_pointer(tensorField, field_fourier, [ res1_red, res(2),res(3),3,3]) ! place a pointer for a complex representation on tensorField
@ -148,7 +151,7 @@ subroutine utilities_init()
call fftw_set_timelimit(fftw_timelimit) ! set timelimit for plan creation call fftw_set_timelimit(fftw_timelimit) ! set timelimit for plan creation
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! creating plans ! creating plans for the convolution
plan_forward = fftw_plan_many_dft_r2c(3, [res(3),res(2) ,res(1)], 9,& ! dimensions, logical length in each dimension in reversed order, no. of transforms plan_forward = fftw_plan_many_dft_r2c(3, [res(3),res(2) ,res(1)], 9,& ! dimensions, logical length in each dimension in reversed order, no. of transforms
field_real, [res(3),res(2) ,res(1)+2_pInt],& ! input data, physical length in each dimension in reversed order field_real, [res(3),res(2) ,res(1)+2_pInt],& ! input data, physical length in each dimension in reversed order
1, res(3)*res(2)*(res(1)+2_pInt),& ! striding, product of physical length in the 3 dimensions 1, res(3)*res(2)*(res(1)+2_pInt),& ! striding, product of physical length in the 3 dimensions
@ -162,7 +165,7 @@ subroutine utilities_init()
1, res(3)*res(2)*(res(1)+2_pInt), fftw_planner_flag) ! striding, product of physical length in the 3 dimensions, planner precision 1, res(3)*res(2)*(res(1)+2_pInt), fftw_planner_flag) ! striding, product of physical length in the 3 dimensions, planner precision
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! depending on (debug) options, allocate more memory and create additional plans ! depending on debug options, allocate more memory and create additional plans
if (debugDivergence) then if (debugDivergence) then
divergence = fftw_alloc_complex(int(res1_red*res(2)*res(3)*3_pInt,C_SIZE_T)) divergence = fftw_alloc_complex(int(res1_red*res(2)*res(3)*3_pInt,C_SIZE_T))
call c_f_pointer(divergence, divergence_real, [ res(1)+2_pInt,res(2),res(3),3]) call c_f_pointer(divergence, divergence_real, [ res(1)+2_pInt,res(2),res(3),3])
@ -267,7 +270,7 @@ end subroutine utilities_updateGamma
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief forward FFT of data in field_real to field_fourier with highest freqs. removed !> @brief forward FFT of data in field_real to field_fourier with highest freqs. removed
!> Does an unweighted FFT transform from real to complex. !> @detailed Does an unweighted FFT transform from real to complex.
!> In case of debugging the FFT, also one component of the tensor (specified by row and column) !> In case of debugging the FFT, also one component of the tensor (specified by row and column)
!> is independetly transformed complex to complex and compared to the whole tensor transform !> is independetly transformed complex to complex and compared to the whole tensor transform
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
@ -326,7 +329,7 @@ end subroutine utilities_FFTforward
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief backward FFT of data in field_fourier to field_real !> @brief backward FFT of data in field_fourier to field_real
!> Does an inverse FFT transform from complex to real !> @detailed Does an inverse FFT transform from complex to real
!> In case of debugging the FFT, also one component of the tensor (specified by row and column) !> In case of debugging the FFT, also one component of the tensor (specified by row and column)
!> is independetly transformed complex to complex and compared to the whole tensor transform !> is independetly transformed complex to complex and compared to the whole tensor transform
!> results is weighted by number of points stored in wgt !> results is weighted by number of points stored in wgt
@ -343,7 +346,7 @@ subroutine utilities_FFTbackward(row,column)
integer(pInt) :: i, j, k, m, n integer(pInt) :: i, j, k, m, n
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! unpack FFT data for conj complex symmetric part ! unpack FFT data for conj complex symmetric part. This data is not transformed when using c2r
if (debugFFTW) then if (debugFFTW) then
scalarField_fourier = field_fourier(1:res1_red,1:res(2),1:res(3),row,column) scalarField_fourier = field_fourier(1:res1_red,1:res(2),1:res(3),row,column)
do i = 0_pInt, res(1)/2_pInt-2_pInt do i = 0_pInt, res(1)/2_pInt-2_pInt
@ -377,6 +380,29 @@ subroutine utilities_FFTbackward(row,column)
field_real = field_real * wgt ! normalize the result by number of elements field_real = field_real * wgt ! normalize the result by number of elements
!--------------------------------------------------------------------------------------------------
! calculate some additional output
! if(debugGeneral) then
! maxCorrectionSkew = 0.0_pReal
! maxCorrectionSym = 0.0_pReal
! temp33_Real = 0.0_pReal
! do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1)
! maxCorrectionSym = max(maxCorrectionSym,&
! maxval(math_symmetric33(field_real(i,j,k,1:3,1:3))))
! maxCorrectionSkew = max(maxCorrectionSkew,&
! maxval(math_skew33(field_real(i,j,k,1:3,1:3))))
! temp33_Real = temp33_Real + field_real(i,j,k,1:3,1:3)
! enddo; enddo; enddo
! write(6,'(a,1x,es11.4)') 'max symmetric correction of deformation =',&
! maxCorrectionSym*wgt
! write(6,'(a,1x,es11.4)') 'max skew correction of deformation =',&
! maxCorrectionSkew*wgt
! write(6,'(a,1x,es11.4)') 'max sym/skew of avg correction = ',&
! maxval(math_symmetric33(temp33_real))/&
! maxval(math_skew33(temp33_real))
! endif
end subroutine utilities_FFTbackward end subroutine utilities_FFTbackward
@ -621,7 +647,7 @@ subroutine utilities_constitutiveResponse(coordinates,F_lastInc,F,temperature,ti
CPFEM_general CPFEM_general
implicit none implicit none
real(pReal), intent(inout), dimension(res(1),res(2),res(3)) :: temperature !< temperature field real(pReal), intent(inout) :: temperature !< temperature (no field)
real(pReal), intent(in), dimension(res(1),res(2),res(3),3) :: coordinates !< coordinates field real(pReal), intent(in), dimension(res(1),res(2),res(3),3) :: coordinates !< coordinates field
real(pReal), intent(in), dimension(3,3,res(1),res(2),res(3)) :: & real(pReal), intent(in), dimension(3,3,res(1),res(2),res(3)) :: &
F_lastInc, & !< target deformation gradient F_lastInc, & !< target deformation gradient
@ -655,7 +681,20 @@ subroutine utilities_constitutiveResponse(coordinates,F_lastInc,F,temperature,ti
calcMode = 2_pInt calcMode = 2_pInt
collectMode = 5_pInt collectMode = 5_pInt
endif endif
!--------------------------------------------------------------------------------------------------
! calculate bounds of det(F) and report
! if(debugGeneral) then
! defgradDetMax = -huge(1.0_pReal)
! defgradDetMin = +huge(1.0_pReal)
! do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1)
! defgradDet = math_det33(F(i,j,k,1:3,1:3))
! defgradDetMax = max(defgradDetMax,defgradDet)
! defgradDetMin = min(defgradDetMin,defgradDet)
! enddo; enddo; enddo
! write(6,'(a,1x,es11.4)') 'max determinant of deformation =', defgradDetMax
! write(6,'(a,1x,es11.4)') 'min determinant of deformation =', defgradDetMin
! endif
if (DebugGeneral) write(6,*) 'collect mode: ', collectMode,' calc mode: ', calcMode if (DebugGeneral) write(6,*) 'collect mode: ', collectMode,' calc mode: ', calcMode
flush(6) flush(6)
@ -664,7 +703,7 @@ subroutine utilities_constitutiveResponse(coordinates,F_lastInc,F,temperature,ti
ielem = ielem + 1_pInt ielem = ielem + 1_pInt
call CPFEM_general(collectMode,& ! collect cycle call CPFEM_general(collectMode,& ! collect cycle
coordinates(i,j,k,1:3), F_lastInc(1:3,1:3,i,j,k),F(1:3,1:3,i,j,k), & coordinates(i,j,k,1:3), F_lastInc(1:3,1:3,i,j,k),F(1:3,1:3,i,j,k), &
temperature(i,j,k),timeinc,ielem,1_pInt,sigma,dsde,P(1:3,1:3,i,j,k),dPdF) temperature,timeinc,ielem,1_pInt,sigma,dsde,P(1:3,1:3,i,j,k),dPdF)
collectMode = 3_pInt collectMode = 3_pInt
enddo; enddo; enddo enddo; enddo; enddo
@ -676,7 +715,7 @@ subroutine utilities_constitutiveResponse(coordinates,F_lastInc,F,temperature,ti
ielem = ielem + 1_pInt ielem = ielem + 1_pInt
call CPFEM_general(calcMode,& ! first element in first iteration retains CPFEM_mode 1, call CPFEM_general(calcMode,& ! first element in first iteration retains CPFEM_mode 1,
coordinates(i,j,k,1:3),F_lastInc(1:3,1:3,i,j,k), F(1:3,1:3,i,j,k), & ! others get 2 (saves winding forward effort) coordinates(i,j,k,1:3),F_lastInc(1:3,1:3,i,j,k), F(1:3,1:3,i,j,k), & ! others get 2 (saves winding forward effort)
temperature(i,j,k),timeinc,ielem,1_pInt,sigma,dsde,P(1:3,1:3,i,j,k),dPdF) temperature,timeinc,ielem,1_pInt,sigma,dsde,P(1:3,1:3,i,j,k),dPdF)
calcMode = 2_pInt calcMode = 2_pInt
C = C + dPdF C = C + dPdF
enddo; enddo; enddo enddo; enddo; enddo