2012-07-19 22:54:56 +05:30
! Copyright 2012 Max-Planck-Institut für Eisenforschung GmbH
!
! This file is part of DAMASK,
! the Düsseldorf Advanced Material Simulation Kit.
!
! DAMASK is free software: you can redistribute it and/or modify
! it under the terms of the GNU General Public License as published by
! the Free Software Foundation, either version 3 of the License, or
! (at your option) any later version.
!
! DAMASK is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! GNU General Public License for more details.
!
! You should have received a copy of the GNU General Public License
! along with DAMASK. If not, see <http://www.gnu.org/licenses/>.
!
!##################################################################################################
!* $Id$
!##################################################################################################
! Material subroutine for BVP solution using spectral method
!
! Run 'DAMASK_spectral.exe --help' to get usage hints
!
! written by P. Eisenlohr,
! F. Roters,
! L. Hantcherli,
! W.A. Counts,
! D.D. Tjahjanto,
! C. Kords,
! M. Diehl,
! R. Lebensohn
!
! MPI fuer Eisenforschung, Duesseldorf
2012-07-20 21:03:13 +05:30
module DAMASK_spectralSolver
2012-07-19 22:54:56 +05:30
use , intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran >4.6 at the moment)
2012-07-20 21:03:13 +05:30
use prec , only : pReal , pInt
use math
use DAMASK_interface , only : &
2012-07-19 22:54:56 +05:30
DAMASK_interface_init , &
loadCaseFile , &
geometryFile , &
getSolverWorkingDirectoryName , &
getSolverJobName , &
appendToOutFile
2012-07-20 21:03:13 +05:30
use debug , only : &
2012-07-19 22:54:56 +05:30
debug_level , &
debug_spectral , &
debug_levelBasic , &
debug_spectralRestart , &
2012-07-20 21:03:13 +05:30
debug_spectralFFTW
use IO
2012-07-23 15:42:31 +05:30
use CPFEM , only : &
CPFEM_general
use numerics , only : &
memory_efficient
2012-07-19 22:54:56 +05:30
implicit none
2012-07-23 15:42:31 +05:30
2012-07-20 21:03:13 +05:30
type solution_t ! mask of stress boundary conditions
logical :: converged = . false .
logical :: regrid = . false .
2012-07-23 15:42:31 +05:30
logical :: term_ill = . false .
2012-07-20 21:03:13 +05:30
end type solution_t
2012-07-23 15:42:31 +05:30
character ( len = 5 ) :: solverType , parameter = 'basic'
real ( pReal ) , dimension ( : , : , : , : , : ) , allocatable :: F , F_lastInc , P
2012-07-20 21:03:13 +05:30
real ( pReal ) , dimension ( : , : , : , : ) , allocatable :: coordinates
real ( pReal ) , dimension ( : , : , : ) , allocatable :: temperature
2012-07-19 22:54:56 +05:30
2012-07-23 15:42:31 +05:30
!--------------------------------------------------------------------------------------------------
! variables storing information for spectral method and FFTW
type ( C_PTR ) :: plan_forward , plan_backward ! plans for fftw
real ( pReal ) , dimension ( 3 , 3 ) :: xiDyad ! product of wave vectors
real ( pReal ) , dimension ( : , : , : , : , : , : , : ) , allocatable :: gamma_hat ! gamma operator (field) for spectral method
real ( pReal ) , dimension ( : , : , : , : ) , allocatable :: xi ! wave vector field for divergence and for gamma operator
real ( pReal ) , dimension ( : , : , : , : , : ) , pointer :: field_real
complex ( pReal ) , dimension ( : , : , : , : , : ) , pointer :: field_fourier
!--------------------------------------------------------------------------------------------------
! debug fftw
type ( C_PTR ) :: plan_scalarField_forth , plan_scalarField_back
complex ( pReal ) , dimension ( : , : , : ) , pointer :: scalarField_real
2012-07-20 21:03:13 +05:30
complex ( pReal ) , dimension ( : , : , : ) , pointer :: scalarField_fourier
2012-07-23 15:42:31 +05:30
!--------------------------------------------------------------------------------------------------
! debug divergence
type ( C_PTR ) :: plan_divergence
2012-07-20 21:03:13 +05:30
real ( pReal ) , dimension ( : , : , : , : ) , pointer :: divergence_real
complex ( pReal ) , dimension ( : , : , : , : ) , pointer :: divergence_fourier
real ( pReal ) , dimension ( : , : , : , : ) , allocatable :: divergence_post
2012-07-23 15:42:31 +05:30
!--------------------------------------------------------------------------------------------------
!variables controlling debugging
logical :: debugGeneral , debugDivergence , debugRestart , debugFFTW
!--------------------------------------------------------------------------------------------------
! variables storing information for spectral method and FFTW
type ( C_PTR ) :: plan_stress , plan_correction ! plans for fftw
!--------------------------------------------------------------------------------------------------
! stress, stiffness and compliance average etc.
real ( pReal ) , dimension ( 3 , 3 ) :: &
F_aim = math_I3 , &
F_aim_lastInc = math_I3
real ( pReal ) , dimension ( 3 , 3 , 3 , 3 ) :: &
C_ref = 0.0_pReal , &
C = 0.0_pReal
real ( pReal ) , dimension ( 3 ) :: geomdim = 0.0_pReal , virt_dim = 0.0_pReal ! physical dimension of volume element per direction
integer ( pInt ) , dimension ( 3 ) :: res = 1_pInt
real ( pReal ) :: wgt
integer ( pInt ) :: res1_red , Npoints
contains
2012-07-20 21:03:13 +05:30
2012-07-23 15:42:31 +05:30
subroutine Solver_Init ( )
2012-07-19 22:54:56 +05:30
2012-07-20 21:03:13 +05:30
use mesh , only : &
mesh_spectral_getResolution , &
mesh_spectral_getDimension
use numerics , only : &
divergence_correction , &
DAMASK_NumThreadsInt , &
fftw_planner_flag , &
fftw_timelimit
use debug , only : &
debug_level , &
debug_spectral , &
debug_levelBasic , &
debug_spectralDivergence , &
debug_spectralRestart , &
2012-07-23 15:42:31 +05:30
debug_spectralFFTW
use FEsolving , only : &
restartInc
2012-07-20 21:03:13 +05:30
implicit none
2012-07-23 15:42:31 +05:30
!--------------------------------------------------------------------------------------------------
2012-07-20 21:03:13 +05:30
! loop variables, convergence etc.
real ( pReal ) , dimension ( 3 , 3 ) , parameter :: ones = 1.0_pReal , zeroes = 0.0_pReal
2012-07-23 15:42:31 +05:30
real ( pReal ) , dimension ( 3 , 3 ) :: temp33_Real
integer ( pInt ) :: i , j , k , l , m , n , q , ielem
real ( pReal ) , dimension ( 3 , 3 , 3 , 3 ) :: dPdF
2012-07-20 21:03:13 +05:30
2012-07-23 15:42:31 +05:30
type ( C_PTR ) :: tensorField ! field in real and fourier space
type ( C_PTR ) :: scalarField_realC , scalarField_fourierC
type ( C_PTR ) :: divergence
2012-07-20 21:03:13 +05:30
2012-07-23 15:42:31 +05:30
integer ( pInt ) , dimension ( 3 ) :: k_s
real ( pReal ) , dimension ( 6 ) :: sigma ! cauchy stress
real ( pReal ) , dimension ( 6 , 6 ) :: dsde
integer ( pInt ) :: ierr
2012-07-19 22:54:56 +05:30
write ( 6 , '(a)' ) ''
2012-07-20 21:03:13 +05:30
write ( 6 , '(a)' ) ' <<<+- DAMASK_spectralSolver init -+>>>'
2012-07-19 22:54:56 +05:30
write ( 6 , '(a)' ) ' $Id$'
#include "compilation_info.f90"
write ( 6 , '(a)' ) ''
!--------------------------------------------------------------------------------------------------
2012-07-23 15:42:31 +05:30
! set debugging parameters
2012-07-19 22:54:56 +05:30
debugGeneral = iand ( debug_level ( debug_spectral ) , debug_levelBasic ) / = 0
debugDivergence = iand ( debug_level ( debug_spectral ) , debug_spectralDivergence ) / = 0
debugRestart = iand ( debug_level ( debug_spectral ) , debug_spectralRestart ) / = 0
debugFFTW = iand ( debug_level ( debug_spectral ) , debug_spectralFFTW ) / = 0
!##################################################################################################
! initialization
!##################################################################################################
2012-07-23 15:42:31 +05:30
res = mesh_spectral_getResolution ( )
geomdim = mesh_spectral_getDimension ( )
res1_red = res ( 1 ) / 2_pInt + 1_pInt
Npoints = res ( 1 ) * res ( 2 ) * res ( 3 )
wgt = 1.0 / real ( Npoints , pReal )
2012-07-19 22:54:56 +05:30
allocate ( F ( res ( 1 ) , res ( 2 ) , res ( 3 ) , 3 , 3 ) , source = 0.0_pReal )
allocate ( F_lastInc ( res ( 1 ) , res ( 2 ) , res ( 3 ) , 3 , 3 ) , source = 0.0_pReal )
2012-07-23 15:42:31 +05:30
allocate ( P ( res ( 1 ) , res ( 2 ) , res ( 3 ) , 3 , 3 ) , source = 0.0_pReal )
2012-07-19 22:54:56 +05:30
allocate ( xi ( 3 , res1_red , res ( 2 ) , res ( 3 ) ) , source = 0.0_pReal )
allocate ( coordinates ( res ( 1 ) , res ( 2 ) , res ( 3 ) , 3 ) , source = 0.0_pReal )
2012-07-23 15:42:31 +05:30
allocate ( temperature ( res ( 1 ) , res ( 2 ) , res ( 3 ) ) , source = 0.0_pReal ) ! start out isothermally
2012-07-19 22:54:56 +05:30
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)
2012-07-23 15:42:31 +05:30
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
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
! general initialization of fftw (see manual on fftw.org for more details)
if ( pReal / = C_DOUBLE . or . pInt / = C_INT ) call IO_error ( error_ID = 808_pInt ) ! check for correct precision in C
!$ if(DAMASK_NumThreadsInt > 0_pInt) then
!$ ierr = fftw_init_threads()
!$ if (ierr == 0_pInt) call IO_error(error_ID = 809_pInt)
!$ call fftw_plan_with_nthreads(DAMASK_NumThreadsInt)
!$ endif
call fftw_set_timelimit ( fftw_timelimit ) ! set timelimit for plan creation
!--------------------------------------------------------------------------------------------------
! creating plans
2012-07-23 15:42:31 +05:30
plan_forward = fftw_plan_many_dft_r2c ( 3 , [ res ( 3 ) , res ( 2 ) , res ( 1 ) ] , 9 , & ! dimensions , 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
2012-07-19 22:54:56 +05:30
1 , res ( 3 ) * res ( 2 ) * ( res ( 1 ) + 2_pInt ) , & ! striding , product of physical lenght in the 3 dimensions
2012-07-23 15:42:31 +05:30
field_fourier , [ res ( 3 ) , res ( 2 ) , res1_red ] , &
2012-07-19 22:54:56 +05:30
1 , res ( 3 ) * res ( 2 ) * res1_red , fftw_planner_flag )
2012-07-23 15:42:31 +05:30
plan_backward = fftw_plan_many_dft_c2r ( 3 , [ res ( 3 ) , res ( 2 ) , res ( 1 ) ] , 9 , &
field_fourier , [ res ( 3 ) , res ( 2 ) , res1_red ] , &
2012-07-19 22:54:56 +05:30
1 , res ( 3 ) * res ( 2 ) * res1_red , &
2012-07-23 15:42:31 +05:30
field_real , [ res ( 3 ) , res ( 2 ) , res ( 1 ) + 2_pInt ] , &
2012-07-19 22:54:56 +05:30
1 , res ( 3 ) * res ( 2 ) * ( res ( 1 ) + 2_pInt ) , fftw_planner_flag )
!--------------------------------------------------------------------------------------------------
! depending on (debug) options, allocate more memory and create additional plans
if ( debugDivergence ) then
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_fourier , [ res1_red , res ( 2 ) , res ( 3 ) , 3 ] )
allocate ( divergence_post ( res ( 1 ) , res ( 2 ) , res ( 3 ) , 3 ) ) ; divergence_post = 0.0_pReal
plan_divergence = fftw_plan_many_dft_c2r ( 3 , [ res ( 3 ) , res ( 2 ) , res ( 1 ) ] , 3 , &
divergence_fourier , [ res ( 3 ) , res ( 2 ) , res1_red ] , &
1 , res ( 3 ) * res ( 2 ) * res1_red , &
divergence_real , [ res ( 3 ) , res ( 2 ) , res ( 1 ) + 2_pInt ] , &
1 , res ( 3 ) * res ( 2 ) * ( res ( 1 ) + 2_pInt ) , fftw_planner_flag )
endif
if ( debugFFTW ) then
scalarField_realC = fftw_alloc_complex ( int ( res ( 1 ) * res ( 2 ) * res ( 3 ) , C_SIZE_T ) ) ! do not do an inplace transform
scalarField_fourierC = fftw_alloc_complex ( int ( res ( 1 ) * res ( 2 ) * res ( 3 ) , C_SIZE_T ) )
call c_f_pointer ( scalarField_realC , scalarField_real , [ res ( 1 ) , res ( 2 ) , res ( 3 ) ] )
call c_f_pointer ( scalarField_fourierC , scalarField_fourier , [ res ( 1 ) , res ( 2 ) , res ( 3 ) ] )
plan_scalarField_forth = fftw_plan_dft_3d ( res ( 3 ) , res ( 2 ) , res ( 1 ) , & !reversed order
scalarField_real , scalarField_fourier , - 1 , fftw_planner_flag )
plan_scalarField_back = fftw_plan_dft_3d ( res ( 3 ) , res ( 2 ) , res ( 1 ) , & !reversed order
scalarField_fourier , scalarField_real , + 1 , fftw_planner_flag )
endif
if ( debugGeneral ) write ( 6 , '(a)' ) 'FFTW initialized'
2012-07-23 15:42:31 +05:30
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
! in case of no restart get reference material stiffness and init fields to no deformation
if ( restartInc == 1_pInt ) then
ielem = 0_pInt
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 ) ; do i = 1_pInt , res ( 1 )
ielem = ielem + 1_pInt
F ( i , j , k , 1 : 3 , 1 : 3 ) = math_I3
F_lastInc ( i , j , k , 1 : 3 , 1 : 3 ) = math_I3
coordinates ( i , j , k , 1 : 3 ) = geomdim / real ( res , pReal ) * real ( [ i , j , k ] , pReal ) - geomdim / real ( 2_pInt * res , pReal )
call CPFEM_general ( 2_pInt , coordinates ( i , j , k , 1 : 3 ) , math_I3 , math_I3 , temperature ( i , j , k ) , &
2012-07-23 15:42:31 +05:30
0.0_pReal , ielem , 1_pInt , sigma , dsde , P ( i , j , k , 1 : 3 , 1 : 3 ) , dPdF )
2012-07-19 22:54:56 +05:30
C = C + dPdF
enddo ; enddo ; enddo
C = C * wgt
C_ref = C
call IO_write_jobBinaryFile ( 777 , 'C_ref' , size ( C_ref ) )
write ( 777 , rec = 1 ) C_ref
close ( 777 )
!--------------------------------------------------------------------------------------------------
! restore deformation gradient and stiffness from saved state
elseif ( restartInc > 1_pInt ) then ! using old values from file
if ( debugRestart ) write ( 6 , '(a,i6,a)' ) 'Reading values of increment ' , &
restartInc - 1_pInt , ' from file'
call IO_read_jobBinaryFile ( 777 , 'convergedSpectralDefgrad' , &
trim ( getSolverJobName ( ) ) , size ( F ) )
read ( 777 , rec = 1 ) F
close ( 777 )
F_lastInc = F
F_aim = 0.0_pReal
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 ) ; do i = 1_pInt , res ( 1 )
F_aim = F_aim + F ( i , j , k , 1 : 3 , 1 : 3 ) ! calculating old average deformation
enddo ; enddo ; enddo
F_aim = F_aim * wgt
F_aim_lastInc = F_aim
coordinates = 0.0 ! change it later!!!
call IO_read_jobBinaryFile ( 777 , 'C_ref' , trim ( getSolverJobName ( ) ) , size ( C_ref ) )
read ( 777 , rec = 1 ) C_ref
close ( 777 )
call IO_read_jobBinaryFile ( 777 , 'C' , trim ( getSolverJobName ( ) ) , size ( C ) )
read ( 777 , rec = 1 ) C
close ( 777 )
endif
!--------------------------------------------------------------------------------------------------
2012-07-23 15:42:31 +05:30
! calculation of discrete angular frequencies, ordered as in FFTW (wrap around)
if ( divergence_correction ) then
do i = 1_pInt , 3_pInt
if ( i / = minloc ( geomdim , 1 ) . and . i / = maxloc ( geomdim , 1 ) ) virt_dim = geomdim / geomdim ( i )
enddo
else
virt_dim = geomdim
endif
do k = 1_pInt , res ( 3 )
k_s ( 3 ) = k - 1_pInt
if ( k > res ( 3 ) / 2_pInt + 1_pInt ) k_s ( 3 ) = k_s ( 3 ) - res ( 3 )
do j = 1_pInt , res ( 2 )
k_s ( 2 ) = j - 1_pInt
if ( j > res ( 2 ) / 2_pInt + 1_pInt ) k_s ( 2 ) = k_s ( 2 ) - res ( 2 )
do i = 1_pInt , res1_red
k_s ( 1 ) = i - 1_pInt
xi ( 1 : 3 , i , j , k ) = real ( k_s , pReal ) / virt_dim
enddo ; enddo ; enddo
!--------------------------------------------------------------------------------------------------
2012-07-19 22:54:56 +05:30
! calculate the gamma operator
if ( memory_efficient ) then ! allocate just single fourth order tensor
allocate ( gamma_hat ( 1 , 1 , 1 , 3 , 3 , 3 , 3 ) , source = 0.0_pReal )
else ! precalculation of gamma_hat field
allocate ( gamma_hat ( res1_red , res ( 2 ) , res ( 3 ) , 3 , 3 , 3 , 3 ) , source = 0.0_pReal )
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 ) ; do i = 1_pInt , res1_red
if ( any ( [ i , j , k ] / = 1_pInt ) ) then ! singular point at xi=(0.0,0.0,0.0) i.e. i=j=k=1
forall ( l = 1_pInt : 3_pInt , m = 1_pInt : 3_pInt ) &
xiDyad ( l , m ) = xi ( l , i , j , k ) * xi ( m , i , j , k )
forall ( l = 1_pInt : 3_pInt , m = 1_pInt : 3_pInt ) &
temp33_Real ( l , m ) = sum ( C_ref ( l , m , 1 : 3 , 1 : 3 ) * xiDyad )
temp33_Real = math_inv33 ( temp33_Real )
2012-07-23 15:42:31 +05:30
forall ( l = 1_pInt : 3_pInt , m = 1_pInt : 3_pInt , n = 1_pInt : 3_pInt , q = 1_pInt : 3_pInt ) &
gamma_hat ( i , j , k , l , m , n , q ) = temp33_Real ( l , n ) * xiDyad ( m , q )
2012-07-19 22:54:56 +05:30
endif
enddo ; enddo ; enddo
gamma_hat ( 1 , 1 , 1 , 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 ) = 0.0_pReal ! singular point at xi=(0.0,0.0,0.0) i.e. i=j=k=1
endif
2012-07-23 15:42:31 +05:30
end subroutine Solver_Init
2012-07-20 21:03:13 +05:30
2012-07-23 15:42:31 +05:30
type ( solution_t ) function solution ( guessmode , timeinc , timeinc_old , P_BC , F_BC , mask_stressVector , velgrad , rotation_BC )
2012-07-20 21:03:13 +05:30
use numerics , only : &
err_div_tol , &
err_stress_tolrel , &
err_stress_tolabs , &
rotation_tol , &
itmax , &
itmin , &
divergence_correction , &
DAMASK_NumThreadsInt , &
fftw_planner_flag , &
fftw_timelimit
2012-07-23 15:42:31 +05:30
use debug , only : &
debug_reset , &
debug_info
!--------------------------------------------------------------------------------------------------
2012-07-20 21:03:13 +05:30
! arrays for mixed boundary conditions
2012-07-23 15:42:31 +05:30
2012-07-20 21:03:13 +05:30
real ( pReal ) , dimension ( 3 , 3 ) , parameter :: ones = 1.0_pReal , zeroes = 0.0_pReal
2012-07-23 15:42:31 +05:30
2012-07-20 21:03:13 +05:30
integer ( pInt ) :: size_reduced = 0_pInt
real ( pReal ) , dimension ( : , : ) , allocatable :: s_reduced , c_reduced ! reduced compliance and stiffness (only for stress BC)
real ( pReal ) , dimension ( 6 ) :: sigma ! cauchy stress
real ( pReal ) , dimension ( 6 , 6 ) :: dsde
real ( pReal ) , dimension ( 9 , 9 ) :: temp99_Real ! compliance and stiffness in matrix notation
2012-07-23 15:42:31 +05:30
real ( pReal ) , dimension ( 3 , 3 , 3 , 3 ) :: &
C_lastInc
integer ( pInt ) :: iter , ielem
2012-07-19 22:54:56 +05:30
2012-07-20 21:03:13 +05:30
!--------------------------------------------------------------------------------------------------
! loop variables, convergence etc.
2012-07-23 15:42:31 +05:30
real ( pReal ) :: timeinc , timeinc_old ! elapsed time, begin of interval, time interval
real ( pReal ) :: guessmode , err_div , err_stress , err_stress_tol
2012-07-20 21:03:13 +05:30
real ( pReal ) , dimension ( 3 , 3 ) :: temp33_Real
2012-07-23 15:42:31 +05:30
integer ( pInt ) :: i , j , k , m , n
integer ( pInt ) :: CPFEM_mode = 1_pInt
logical :: errmatinv , restartWrite
2012-07-20 21:03:13 +05:30
real ( pReal ) :: defgradDet
2012-07-19 22:54:56 +05:30
2012-07-23 15:42:31 +05:30
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
2012-07-20 21:03:13 +05:30
! variables for additional output of divergence calculations
2012-07-23 15:42:31 +05:30
2012-07-20 21:03:13 +05:30
real ( pReal ) , dimension ( 3 , 3 ) :: &
P_av , &
mask_stress , &
mask_defgrad , &
deltaF_aim , &
F_aim_lab , &
F_aim_lab_lastIter , &
P_av_lab
2012-07-23 15:42:31 +05:30
real ( pReal ) , dimension ( 3 , 3 , 3 , 3 ) :: &
2012-07-20 21:03:13 +05:30
dPdF , &
C = 0.0_pReal , &
2012-07-23 15:42:31 +05:30
S_lastInc
2012-07-20 21:03:13 +05:30
2012-07-23 15:42:31 +05:30
logical :: velgrad
real ( pReal ) , dimension ( 3 , 3 ) :: P_BC , F_BC , rotation_BC
logical , dimension ( 9 ) :: mask_stressVector
real ( pReal ) :: defgradDetMax , defgradDetMin
mask_stress = merge ( ones , zeroes , reshape ( mask_stressVector , [ 3 , 3 ] ) )
mask_defgrad = merge ( zeroes , ones , reshape ( mask_stressVector , [ 3 , 3 ] ) )
size_reduced = int ( count ( mask_stressVector ) , pInt )
2012-07-20 21:03:13 +05:30
allocate ( c_reduced ( size_reduced , size_reduced ) , source = 0.0_pReal )
allocate ( s_reduced ( size_reduced , size_reduced ) , source = 0.0_pReal )
2012-07-23 15:42:31 +05:30
if ( velgrad ) then ! calculate deltaF_aim from given L and current F
deltaF_aim = timeinc * mask_defgrad * math_mul33x33 ( F_BC , F_aim )
else ! deltaF_aim = fDot *timeinc where applicable
deltaF_aim = timeinc * mask_defgrad * F_BC
endif
!--------------------------------------------------------------------------------------------------
! winding forward of deformation aim in loadcase system
temp33_Real = F_aim
F_aim = F_aim &
+ guessmode * mask_stress * ( F_aim - F_aim_lastInc ) * timeinc / timeinc_old &
+ deltaF_aim
F_aim_lastInc = temp33_Real
!--------------------------------------------------------------------------------------------------
! winding forward of deformation aim in loadcase system
temp33_Real = F_aim
F_aim = F_aim &
+ guessmode * mask_stress * ( F_aim - F_aim_lastInc ) * timeinc / timeinc_old &
+ deltaF_aim
F_aim_lastInc = temp33_Real
!--------------------------------------------------------------------------------------------------
! update local deformation gradient and coordinates
deltaF_aim = math_rotate_backward33 ( deltaF_aim , rotation_BC )
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 ) ; do i = 1_pInt , res ( 1 )
temp33_Real = F ( i , j , k , 1 : 3 , 1 : 3 )
F ( i , j , k , 1 : 3 , 1 : 3 ) = F ( i , j , k , 1 : 3 , 1 : 3 ) & ! decide if guessing along former trajectory or apply homogeneous addon
+ guessmode * ( F ( i , j , k , 1 : 3 , 1 : 3 ) - F_lastInc ( i , j , k , 1 : 3 , 1 : 3 ) ) & ! guessing...
* timeinc / timeinc_old &
+ ( 1.0_pReal - guessmode ) * deltaF_aim ! if not guessing, use prescribed average deformation where applicable
F_lastInc ( i , j , k , 1 : 3 , 1 : 3 ) = temp33_Real
enddo ; enddo ; enddo
call deformed_fft ( res , geomdim , math_rotate_backward33 ( F_aim , rotation_BC ) , & ! calculate current coordinates
1.0_pReal , F_lastInc , coordinates )
2012-07-20 21:03:13 +05:30
2012-07-23 15:42:31 +05:30
2012-07-20 21:03:13 +05:30
!--------------------------------------------------------------------------------------------------
2012-07-19 22:54:56 +05:30
! calculate reduced compliance
if ( size_reduced > 0_pInt ) then ! calculate compliance in case stress BC is applied
2012-07-23 15:42:31 +05:30
C_lastInc = math_rotate_forward3333 ( C , rotation_BC ) ! calculate stiffness from former inc
2012-07-19 22:54:56 +05:30
temp99_Real = math_Plain3333to99 ( C_lastInc )
k = 0_pInt ! build reduced stiffness
do n = 1_pInt , 9_pInt
2012-07-23 15:42:31 +05:30
if ( mask_stressVector ( n ) ) then
2012-07-19 22:54:56 +05:30
k = k + 1_pInt
j = 0_pInt
do m = 1_pInt , 9_pInt
2012-07-23 15:42:31 +05:30
if ( mask_stressVector ( m ) ) then
2012-07-19 22:54:56 +05:30
j = j + 1_pInt
c_reduced ( k , j ) = temp99_Real ( n , m )
endif ; enddo ; endif ; enddo
call math_invert ( size_reduced , c_reduced , s_reduced , i , errmatinv ) ! invert reduced stiffness
if ( errmatinv ) call IO_error ( error_ID = 400_pInt )
temp99_Real = 0.0_pReal ! build full compliance
k = 0_pInt
do n = 1_pInt , 9_pInt
2012-07-23 15:42:31 +05:30
if ( mask_stressVector ( n ) ) then
2012-07-19 22:54:56 +05:30
k = k + 1_pInt
j = 0_pInt
do m = 1_pInt , 9_pInt
2012-07-23 15:42:31 +05:30
if ( mask_stressVector ( m ) ) then
2012-07-19 22:54:56 +05:30
j = j + 1_pInt
temp99_Real ( n , m ) = s_reduced ( k , j )
endif ; enddo ; endif ; enddo
S_lastInc = ( math_Plain99to3333 ( temp99_Real ) )
endif
2012-07-23 15:42:31 +05:30
2012-07-19 22:54:56 +05:30
iter = 0_pInt
err_div = huge ( err_div_tol ) ! go into loop
2012-07-23 15:42:31 +05:30
2012-07-20 21:03:13 +05:30
!##################################################################################################
2012-07-19 22:54:56 +05:30
! convergence loop (looping over iterations)
!##################################################################################################
2012-07-20 21:03:13 +05:30
2012-07-19 22:54:56 +05:30
do while ( ( iter < itmax . and . ( err_div > err_div_tol . or . err_stress > err_stress_tol ) ) &
. or . iter < itmin )
iter = iter + 1_pInt
2012-07-23 15:42:31 +05:30
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
P_av = math_rotate_forward33 ( P_av_lab , rotation_BC )
write ( 6 , '(a,/,3(3(f12.7,1x)/))' , advance = 'no' ) 'Piola-Kirchhoff stress / MPa =' , &
math_transpose33 ( P_av ) / 1.e6_pReal
!--------------------------------------------------------------------------------------------------
! stress BC handling
if ( size_reduced > 0_pInt ) then ! calculate stress BC if applied
err_stress = maxval ( abs ( mask_stress * ( P_av - P_BC ) ) ) ! maximum deviaton (tensor norm not applicable)
err_stress_tol = min ( maxval ( abs ( P_av ) ) * err_stress_tolrel , err_stress_tolabs ) ! don't use any tensor norm for the relative criterion because the comparison should be coherent
write ( 6 , '(a)' ) ''
write ( 6 , '(a)' ) '... correcting deformation gradient to fulfill BCs ...............'
write ( 6 , '(a,f6.2,a,es11.4,a)' ) 'error stress = ' , err_stress / err_stress_tol , &
' (' , err_stress , ' Pa)'
F_aim = F_aim - math_mul3333xx33 ( S_lastInc , ( ( P_av - P_BC ) ) ) ! residual on given stress components
write ( 6 , '(a,1x,es11.4)' ) 'determinant of new deformation = ' , math_det33 ( F_aim )
else
err_stress_tol = + huge ( 1.0_pReal )
endif
F_aim_lab = math_rotate_backward33 ( F_aim , rotation_BC ) ! boundary conditions from load frame into lab (Fourier) frame
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
! report begin of new iteration
write ( 6 , '(a)' ) ''
write ( 6 , '(a)' ) '=================================================================='
2012-07-23 15:42:31 +05:30
write ( 6 , '(3(a,i6.6))' ) ' @ Iter. ' , itmin , ' < ' , iter , ' < ' , itmax
2012-07-19 22:54:56 +05:30
write ( 6 , '(a,/,3(3(f12.7,1x)/))' , advance = 'no' ) 'deformation gradient aim =' , &
math_transpose33 ( F_aim )
write ( 6 , '(a)' ) ''
write ( 6 , '(a)' ) '... update stress field P(F) .....................................'
if ( restartWrite ) write ( 6 , '(a)' ) 'writing restart info for last increment'
2012-07-23 15:42:31 +05:30
if ( restartWrite ) then
write ( 6 , '(a)' ) 'writing converged results for restart'
call IO_write_jobBinaryFile ( 777 , 'convergedSpectralDefgrad' , size ( F_lastInc ) ) ! writing deformation gradient field to file
write ( 777 , rec = 1 ) F_LastInc
close ( 777 )
call IO_write_jobBinaryFile ( 777 , 'C' , size ( C ) )
write ( 777 , rec = 1 ) C
close ( 777 )
endif
F_aim_lab_lastIter = math_rotate_backward33 ( F_aim , rotation_BC )
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
! evaluate constitutive response
ielem = 0_pInt
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 ) ; do i = 1_pInt , res ( 1 )
ielem = ielem + 1_pInt
call CPFEM_general ( 3_pInt , & ! collect cycle
coordinates ( i , j , k , 1 : 3 ) , F_lastInc ( i , j , k , 1 : 3 , 1 : 3 ) , F ( i , j , k , 1 : 3 , 1 : 3 ) , &
temperature ( i , j , k ) , timeinc , ielem , 1_pInt , sigma , dsde , &
2012-07-23 15:42:31 +05:30
P ( i , j , k , 1 : 3 , 1 : 3 ) , dPdF )
2012-07-19 22:54:56 +05:30
enddo ; enddo ; enddo
2012-07-23 15:42:31 +05:30
P = 0.0_pReal ! needed because of the padding for FFTW
2012-07-19 22:54:56 +05:30
C = 0.0_pReal
2012-07-23 15:42:31 +05:30
P_av_lab = 0.0_pReal
2012-07-19 22:54:56 +05:30
ielem = 0_pInt
call debug_reset ( )
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 ) ; do i = 1_pInt , res ( 1 )
ielem = ielem + 1_pInt
call CPFEM_general ( CPFEM_mode , & ! first element in first iteration retains CPFEM_mode 1,
coordinates ( i , j , k , 1 : 3 ) , F_lastInc ( i , j , k , 1 : 3 , 1 : 3 ) , F ( i , j , k , 1 : 3 , 1 : 3 ) , & ! others get 2 (saves winding forward effort)
temperature ( i , j , k ) , timeinc , ielem , 1_pInt , sigma , dsde , &
2012-07-23 15:42:31 +05:30
P ( i , j , k , 1 : 3 , 1 : 3 ) , dPdF )
2012-07-19 22:54:56 +05:30
CPFEM_mode = 2_pInt
C = C + dPdF
2012-07-23 15:42:31 +05:30
P_av_lab = P_av_lab + P ( i , j , k , 1 : 3 , 1 : 3 )
2012-07-19 22:54:56 +05:30
enddo ; enddo ; enddo
call debug_info ( )
2012-07-23 15:42:31 +05:30
restartWrite = . false .
P_av_lab = P_av_lab * wgt
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
P_av = math_rotate_forward33 ( P_av_lab , rotation_BC )
write ( 6 , '(a,/,3(3(f12.7,1x)/))' , advance = 'no' ) 'Piola-Kirchhoff stress / MPa =' , &
math_transpose33 ( P_av ) / 1.e6_pReal
!--------------------------------------------------------------------------------------------------
! stress BC handling
if ( size_reduced > 0_pInt ) then ! calculate stress BC if applied
err_stress = maxval ( abs ( mask_stress * ( P_av - P_BC ) ) ) ! maximum deviaton (tensor norm not applicable)
err_stress_tol = min ( maxval ( abs ( P_av ) ) * err_stress_tolrel , err_stress_tolabs ) ! don't use any tensor norm for the relative criterion because the comparison should be coherent
write ( 6 , '(a)' ) ''
write ( 6 , '(a)' ) '... correcting deformation gradient to fulfill BCs ...............'
write ( 6 , '(a,f6.2,a,es11.4,a)' ) 'error stress = ' , err_stress / err_stress_tol , &
' (' , err_stress , ' Pa)'
F_aim = F_aim - math_mul3333xx33 ( S_lastInc , ( ( P_av - P_BC ) ) ) ! residual on given stress components
write ( 6 , '(a,1x,es11.4)' ) 'determinant of new deformation = ' , math_det33 ( F_aim )
else
err_stress_tol = + huge ( 1.0_pReal )
endif
F_aim_lab = math_rotate_backward33 ( F_aim , rotation_BC )
field_real ( 1 : res ( 1 ) , 1 : res ( 2 ) , 1 : res ( 3 ) , 1 : 3 , 1 : 3 ) = P
call convolution ( . True . , F_aim_lab_lastIter - F_aim_lab )
!--------------------------------------------------------------------------------------------------
! updated deformation gradient
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 ) ; do i = 1_pInt , res ( 1 )
F ( i , j , k , 1 : 3 , 1 : 3 ) = F ( i , j , k , 1 : 3 , 1 : 3 ) - field_real ( i , j , k , 1 : 3 , 1 : 3 ) * wgt ! F(x)^(n+1) = F(x)^(n) + correction; *wgt: correcting for missing normalization
enddo ; enddo ; enddo
!--------------------------------------------------------------------------------------------------
! 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
enddo ! end looping when convergency is achieved
CPFEM_mode = 1_pInt ! winding forward
C = C * wgt
deallocate ( c_reduced )
deallocate ( s_reduced )
end function solution
subroutine convolution ( calcDivergence , field_aim )
real ( pReal ) :: err_div
real ( pReal ) , dimension ( 3 , 3 ) :: temp33_Real
integer ( pInt ) :: i , j , k , l , m , n , q
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
2012-07-23 15:42:31 +05:30
!variables for additional output due to general debugging
real ( pReal ) :: maxCorrectionSym , maxCorrectionSkew
logical :: calcDivergence
real ( pReal ) , dimension ( 3 , 3 ) :: field_avg , field_aim
integer ( pInt ) :: row , column
real ( pReal ) :: field_av_L2 , err_div_RMS , err_real_div_RMS , err_post_div_RMS , &
err_div_max , err_real_div_max
complex ( pReal ) , dimension ( 3 ) :: temp3_complex
complex ( pReal ) , dimension ( 3 , 3 ) :: temp33_complex
!--------------------------------------------------------------------------------------------------
! actual spectral method
write ( 6 , '(a)' ) ''
write ( 6 , '(a)' ) '... doing convolution .................'
!--------------------------------------------------------------------------------------------------
2012-07-19 22:54:56 +05:30
! copy one component of the stress field to to a single FT and check for mismatch
if ( debugFFTW ) then
2012-07-23 15:42:31 +05:30
row = 3 ! (mod(totalIncsCounter+iter-2_pInt,9_pInt))/3_pInt + 1_pInt ! go through the elements of the tensors, controlled by totalIncsCounter and iter, starting at 1
column = 3 !(mod(totalIncsCounter+iter-2_pInt,3_pInt)) + 1_pInt
2012-07-19 22:54:56 +05:30
scalarField_real ( 1 : res ( 1 ) , 1 : res ( 2 ) , 1 : res ( 3 ) ) = & ! store the selected component
2012-07-23 15:42:31 +05:30
cmplx ( field_real ( 1 : res ( 1 ) , 1 : res ( 2 ) , 1 : res ( 3 ) , row , column ) , 0.0_pReal , pReal )
2012-07-19 22:54:56 +05:30
endif
!--------------------------------------------------------------------------------------------------
! call function to calculate divergence from math (for post processing) to check results
if ( debugDivergence ) &
call divergence_fft ( res , virt_dim , 3_pInt , &
2012-07-23 15:42:31 +05:30
field_real ( 1 : res ( 1 ) , 1 : res ( 2 ) , 1 : res ( 3 ) , 1 : 3 , 1 : 3 ) , divergence_post ) ! padding
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
! doing the FT because it simplifies calculation of average stress in real space also
2012-07-23 15:42:31 +05:30
call fftw_execute_dft_r2c ( plan_forward , field_real , field_fourier )
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
! comparing 1 and 3x3 FT results
if ( debugFFTW ) then
call fftw_execute_dft ( plan_scalarField_forth , scalarField_real , scalarField_fourier )
write ( 6 , '(a,i1,1x,i1)' ) 'checking FT results of compontent ' , row , column
write ( 6 , '(a,2(es11.4,1x))' ) 'max FT relative error = ' , &
maxval ( real ( ( scalarField_fourier ( 1 : res1_red , 1 : res ( 2 ) , 1 : res ( 3 ) ) - &
2012-07-23 15:42:31 +05:30
field_fourier ( 1 : res1_red , 1 : res ( 2 ) , 1 : res ( 3 ) , row , column ) ) / &
2012-07-19 22:54:56 +05:30
scalarField_fourier ( 1 : res1_red , 1 : res ( 2 ) , 1 : res ( 3 ) ) ) ) , &
maxval ( aimag ( ( scalarField_fourier ( 1 : res1_red , 1 : res ( 2 ) , 1 : res ( 3 ) ) - &
2012-07-23 15:42:31 +05:30
field_fourier ( 1 : res1_red , 1 : res ( 2 ) , 1 : res ( 3 ) , row , column ) ) / &
2012-07-19 22:54:56 +05:30
scalarField_fourier ( 1 : res1_red , 1 : res ( 2 ) , 1 : res ( 3 ) ) ) )
endif
!--------------------------------------------------------------------------------------------------
! removing highest frequencies
2012-07-23 15:42:31 +05:30
field_fourier ( res1_red , 1 : res ( 2 ) , 1 : res ( 3 ) , 1 : 3 , 1 : 3 ) &
2012-07-19 22:54:56 +05:30
= cmplx ( 0.0_pReal , 0.0_pReal , pReal )
2012-07-23 15:42:31 +05:30
field_fourier ( 1 : res1_red , res ( 2 ) / 2_pInt + 1_pInt , 1 : res ( 3 ) , 1 : 3 , 1 : 3 ) &
2012-07-19 22:54:56 +05:30
= cmplx ( 0.0_pReal , 0.0_pReal , pReal )
if ( res ( 3 ) > 1_pInt ) &
2012-07-23 15:42:31 +05:30
field_fourier ( 1 : res1_red , 1 : res ( 2 ) , res ( 3 ) / 2_pInt + 1_pInt , 1 : 3 , 1 : 3 ) &
2012-07-19 22:54:56 +05:30
= cmplx ( 0.0_pReal , 0.0_pReal , pReal )
!--------------------------------------------------------------------------------------------------
! calculating RMS divergence criterion in Fourier space
2012-07-23 15:42:31 +05:30
if ( calcDivergence ) then
field_avg = real ( field_fourier ( 1 , 1 , 1 , 1 : 3 , 1 : 3 ) , pReal ) * wgt
field_av_L2 = sqrt ( maxval ( math_eigenvalues33 ( math_mul33x33 ( field_avg , & ! L_2 norm of average stress (http://mathworld.wolfram.com/SpectralNorm.html)
math_transpose33 ( field_avg ) ) ) ) )
2012-07-19 22:54:56 +05:30
err_div_RMS = 0.0_pReal
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 )
do i = 2_pInt , res1_red - 1_pInt ! Has somewhere a conj. complex counterpart. Therefore count it twice.
err_div_RMS = err_div_RMS &
2012-07-23 15:42:31 +05:30
+ 2.0_pReal * ( sum ( real ( math_mul33x3_complex ( field_fourier ( i , j , k , 1 : 3 , 1 : 3 ) , & ! (sqrt(real(a)**2 + aimag(a)**2))**2 = real(a)**2 + aimag(a)**2. do not take square root and square again
2012-07-19 22:54:56 +05:30
xi ( 1 : 3 , i , j , k ) ) * TWOPIIMG ) ** 2.0_pReal ) & ! --> sum squared L_2 norm of vector
2012-07-23 15:42:31 +05:30
+ sum ( aimag ( math_mul33x3_complex ( field_fourier ( i , j , k , 1 : 3 , 1 : 3 ) , &
2012-07-19 22:54:56 +05:30
xi ( 1 : 3 , i , j , k ) ) * TWOPIIMG ) ** 2.0_pReal ) )
enddo
err_div_RMS = err_div_RMS & ! Those two layers (DC and Nyquist) do not have a conjugate complex counterpart
2012-07-23 15:42:31 +05:30
+ sum ( real ( math_mul33x3_complex ( field_fourier ( 1 , j , k , 1 : 3 , 1 : 3 ) , &
2012-07-19 22:54:56 +05:30
xi ( 1 : 3 , 1 , j , k ) ) * TWOPIIMG ) ** 2.0_pReal ) &
2012-07-23 15:42:31 +05:30
+ sum ( aimag ( math_mul33x3_complex ( field_fourier ( 1 , j , k , 1 : 3 , 1 : 3 ) , &
2012-07-19 22:54:56 +05:30
xi ( 1 : 3 , 1 , j , k ) ) * TWOPIIMG ) ** 2.0_pReal ) &
2012-07-23 15:42:31 +05:30
+ sum ( real ( math_mul33x3_complex ( field_fourier ( res1_red , j , k , 1 : 3 , 1 : 3 ) , &
2012-07-19 22:54:56 +05:30
xi ( 1 : 3 , res1_red , j , k ) ) * TWOPIIMG ) ** 2.0_pReal ) &
2012-07-23 15:42:31 +05:30
+ sum ( aimag ( math_mul33x3_complex ( field_fourier ( res1_red , j , k , 1 : 3 , 1 : 3 ) , &
2012-07-19 22:54:56 +05:30
xi ( 1 : 3 , res1_red , j , k ) ) * TWOPIIMG ) ** 2.0_pReal )
enddo ; enddo
err_div_RMS = sqrt ( err_div_RMS ) * wgt ! RMS in real space calculated with Parsevals theorem from Fourier space
2012-07-23 15:42:31 +05:30
err_div = err_div_RMS / field_av_L2 ! criterion to stop iterations
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
! calculate additional divergence criteria and report
if ( debugDivergence ) then ! calculate divergence again
err_div_max = 0.0_pReal
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 ) ; do i = 1_pInt , res1_red
2012-07-23 15:42:31 +05:30
temp3_Complex = math_mul33x3_complex ( field_fourier ( i , j , k , 1 : 3 , 1 : 3 ) * wgt , & ! weighting P_fourier
2012-07-19 22:54:56 +05:30
xi ( 1 : 3 , i , j , k ) ) * TWOPIIMG
err_div_max = max ( err_div_max , sum ( abs ( temp3_Complex ) ** 2.0_pReal ) )
divergence_fourier ( i , j , k , 1 : 3 ) = temp3_Complex ! need divergence NOT squared
enddo ; enddo ; enddo
call fftw_execute_dft_c2r ( plan_divergence , divergence_fourier , divergence_real ) ! already weighted
err_real_div_RMS = 0.0_pReal
err_post_div_RMS = 0.0_pReal
err_real_div_max = 0.0_pReal
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 ) ; do i = 1_pInt , res ( 1 )
err_real_div_RMS = err_real_div_RMS + sum ( divergence_real ( i , j , k , 1 : 3 ) ** 2.0_pReal ) ! avg of squared L_2 norm of div(stress) in real space
err_post_div_RMS = err_post_div_RMS + sum ( divergence_post ( i , j , k , 1 : 3 ) ** 2.0_pReal ) ! avg of squared L_2 norm of div(stress) in real space
err_real_div_max = max ( err_real_div_max , sum ( divergence_real ( i , j , k , 1 : 3 ) ** 2.0_pReal ) ) ! max of squared L_2 norm of div(stress) in real space
enddo ; enddo ; enddo
err_real_div_RMS = sqrt ( wgt * err_real_div_RMS ) ! RMS in real space
err_post_div_RMS = sqrt ( wgt * err_post_div_RMS ) ! RMS in real space
err_real_div_max = sqrt ( err_real_div_max ) ! max in real space
err_div_max = sqrt ( err_div_max ) ! max in Fourier space
write ( 6 , '(a,es11.4)' ) 'error divergence FT RMS = ' , err_div_RMS
write ( 6 , '(a,es11.4)' ) 'error divergence Real RMS = ' , err_real_div_RMS
write ( 6 , '(a,es11.4)' ) 'error divergence post RMS = ' , err_post_div_RMS
write ( 6 , '(a,es11.4)' ) 'error divergence FT max = ' , err_div_max
write ( 6 , '(a,es11.4)' ) 'error divergence Real max = ' , err_real_div_max
endif
2012-07-23 15:42:31 +05:30
! write(6,'(a,f6.2,a,es11.4,a)') 'error divergence = ', err_div/err_div_tol,&
! ' (',err_div,' N/m³)'
end if
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
! to the actual spectral method calculation (mechanical equilibrium)
if ( memory_efficient ) then ! memory saving version, on-the-fly calculation of gamma_hat
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 ) ; do i = 1_pInt , res1_red
if ( any ( [ i , j , k ] / = 1_pInt ) ) then ! singular point at xi=(0.0,0.0,0.0) i.e. i=j=k=1
forall ( l = 1_pInt : 3_pInt , m = 1_pInt : 3_pInt ) &
xiDyad ( l , m ) = xi ( l , i , j , k ) * xi ( m , i , j , k )
forall ( l = 1_pInt : 3_pInt , m = 1_pInt : 3_pInt ) &
temp33_Real ( l , m ) = sum ( C_ref ( l , m , 1 : 3 , 1 : 3 ) * xiDyad )
temp33_Real = math_inv33 ( temp33_Real )
2012-07-23 15:42:31 +05:30
forall ( l = 1_pInt : 3_pInt , m = 1_pInt : 3_pInt , n = 1_pInt : 3_pInt , q = 1_pInt : 3_pInt ) &
gamma_hat ( 1 , 1 , 1 , l , m , n , q ) = temp33_Real ( l , n ) * xiDyad ( m , q )
2012-07-19 22:54:56 +05:30
forall ( l = 1_pInt : 3_pInt , m = 1_pInt : 3_pInt ) &
temp33_Complex ( l , m ) = sum ( gamma_hat ( 1 , 1 , 1 , l , m , 1 : 3 , 1 : 3 ) * &
2012-07-23 15:42:31 +05:30
field_fourier ( i , j , k , 1 : 3 , 1 : 3 ) )
field_fourier ( i , j , k , 1 : 3 , 1 : 3 ) = temp33_Complex
2012-07-19 22:54:56 +05:30
endif
enddo ; enddo ; enddo
else ! use precalculated gamma-operator
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 ) ; do i = 1_pInt , res1_red
forall ( m = 1_pInt : 3_pInt , n = 1_pInt : 3_pInt ) &
temp33_Complex ( m , n ) = sum ( gamma_hat ( i , j , k , m , n , 1 : 3 , 1 : 3 ) * &
2012-07-23 15:42:31 +05:30
field_fourier ( i , j , k , 1 : 3 , 1 : 3 ) )
field_fourier ( i , j , k , 1 : 3 , 1 : 3 ) = temp33_Complex
2012-07-19 22:54:56 +05:30
enddo ; enddo ; enddo
endif
2012-07-23 15:42:31 +05:30
field_fourier ( 1 , 1 , 1 , 1 : 3 , 1 : 3 ) = cmplx ( field_aim , 0.0_pReal , pReal ) ! singular point at xi=(0.0,0.0,0.0) i.e. i=j=k=1
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
! comparing 1 and 3x3 inverse FT results
if ( debugFFTW ) then
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 ) ; do i = 1_pInt , res1_red
2012-07-23 15:42:31 +05:30
scalarField_fourier ( i , j , k ) = field_fourier ( i , j , k , row , column )
2012-07-19 22:54:56 +05:30
enddo ; enddo ; enddo
do i = 0_pInt , res ( 1 ) / 2_pInt - 2_pInt ! unpack fft data for conj complex symmetric part
m = 1_pInt
do k = 1_pInt , res ( 3 )
n = 1_pInt
do j = 1_pInt , res ( 2 )
scalarField_fourier ( res ( 1 ) - i , j , k ) = conjg ( scalarField_fourier ( 2 + i , n , m ) )
if ( n == 1_pInt ) n = res ( 2 ) + 1_pInt
n = n - 1_pInt
enddo
if ( m == 1_pInt ) m = res ( 3 ) + 1_pInt
m = m - 1_pInt
enddo ; enddo
endif
!--------------------------------------------------------------------------------------------------
! doing the inverse FT
2012-07-23 15:42:31 +05:30
call fftw_execute_dft_c2r ( plan_backward , field_fourier , field_real ) ! back transform of fluct deformation gradient
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
! comparing 1 and 3x3 inverse FT results
if ( debugFFTW ) then
write ( 6 , '(a,i1,1x,i1)' ) 'checking iFT results of compontent ' , row , column
call fftw_execute_dft ( plan_scalarField_back , scalarField_fourier , scalarField_real )
write ( 6 , '(a,es11.4)' ) 'max iFT relative error = ' , &
maxval ( ( real ( scalarField_real ( 1 : res ( 1 ) , 1 : res ( 2 ) , 1 : res ( 3 ) ) ) - &
2012-07-23 15:42:31 +05:30
field_real ( 1 : res ( 1 ) , 1 : res ( 2 ) , 1 : res ( 3 ) , row , column ) ) / &
2012-07-19 22:54:56 +05:30
real ( scalarField_real ( 1 : res ( 1 ) , 1 : res ( 2 ) , 1 : res ( 3 ) ) ) )
endif
!--------------------------------------------------------------------------------------------------
! 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 , &
2012-07-23 15:42:31 +05:30
maxval ( math_symmetric33 ( field_real ( i , j , k , 1 : 3 , 1 : 3 ) ) ) )
2012-07-19 22:54:56 +05:30
maxCorrectionSkew = max ( maxCorrectionSkew , &
2012-07-23 15:42:31 +05:30
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 )
2012-07-19 22:54:56 +05:30
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
2012-07-23 15:42:31 +05:30
end subroutine convolution
2012-07-20 21:03:13 +05:30
end module DAMASK_spectralSolver