2012-03-07 23:07:40 +05:30
! Copyright 2011 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/>.
!
!##################################################################################################
2012-03-30 01:24:31 +05:30
!* $Id$
2012-03-07 23:07:40 +05:30
!##################################################################################################
! 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-04-11 22:58:08 +05:30
#include "spectral_quit.f90"
2012-03-07 23:07:40 +05:30
2012-04-11 22:58:08 +05:30
program DAMASK_spectral_AL
2012-03-07 23:07:40 +05:30
use , intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran >4.6 at the moment)
2012-04-11 22:58:08 +05:30
use DAMASK_interface , only : &
DAMASK_interface_init , &
getLoadcaseName , &
getSolverWorkingDirectoryName , &
getSolverJobName , &
getModelName , &
inputFileExtension
use prec , only : &
pInt , &
pReal , &
DAMASK_NaN
use IO , only : &
IO_isBlank , &
IO_open_file , &
IO_stringPos , &
IO_stringValue , &
IO_floatValue , &
IO_intValue , &
IO_error , &
IO_lc , &
IO_read_jobBinaryFile , &
IO_write_jobBinaryFile
use debug , only : &
debug_what , &
debug_spectral , &
debug_levelBasic , &
debug_spectralDivergence , &
debug_spectralRestart , &
debug_spectralFFTW , &
debug_reset , &
debug_info
2012-03-07 23:07:40 +05:30
use math
2012-04-11 22:58:08 +05:30
use mesh , only : &
mesh_spectral_getResolution , &
mesh_spectral_getDimension , &
mesh_spectral_getHomogenization
use CPFEM , only : &
CPFEM_general , &
CPFEM_initAll
use FEsolving , only : &
restartWrite , &
restartInc
use numerics , only : &
err_stress_tolrel , &
err_stress_tolabs , &
rotation_tol , &
itmax , &
itmin , &
memory_efficient , &
update_gamma , &
divergence_correction , &
DAMASK_NumThreadsInt , &
fftw_planner_flag , &
fftw_timelimit
use homogenization , only : &
materialpoint_sizeResults , &
materialpoint_results
implicit none
2012-03-07 23:07:40 +05:30
!--------------------------------------------------------------------------------------------------
! variables to read from load case and geom file
real ( pReal ) , dimension ( 9 ) :: temp_valueVector ! stores information temporarily from loadcase file
logical , dimension ( 9 ) :: temp_maskVector
integer ( pInt ) , parameter :: maxNchunksLoadcase = ( 1_pInt + 9_pInt ) * 3_pInt + & ! deformation, rotation, and stress
( 1_pInt + 1_pInt ) * 5_pInt + & ! time, (log)incs, temp, restartfrequency, and outputfrequency
1_pInt , & ! dropguessing
maxNchunksGeom = 7_pInt , & ! 4 identifiers, 3 values
myUnit = 234_pInt
integer ( pInt ) , dimension ( 1_pInt + maxNchunksLoadcase * 2_pInt ) :: positions ! this is longer than needed for geometry parsing
2012-04-11 22:58:08 +05:30
integer ( pInt ) :: &
N_l = 0_pInt , &
N_t = 0_pInt , &
N_n = 0_pInt , &
N_Fdot = 0_pInt
character ( len = 1024 ) :: line
2012-03-07 23:07:40 +05:30
!--------------------------------------------------------------------------------------------------
! variable storing information from load case file
type bc_type
real ( pReal ) , dimension ( 3 , 3 ) :: deformation = 0.0_pReal , & ! applied velocity gradient or time derivative of deformation gradient
2012-03-20 23:31:31 +05:30
P = 0.0_pReal , & ! stress BC (if applicable)
2012-03-07 23:07:40 +05:30
rotation = math_I3 ! rotation of BC (if applicable)
real ( pReal ) :: time = 0.0_pReal , & ! length of increment
temperature = 30 0.0_pReal ! isothermal starting conditions
integer ( pInt ) :: incs = 0_pInt , & ! number of increments
outputfrequency = 1_pInt , & ! frequency of result writes
restartfrequency = 0_pInt , & ! frequency of restart writes
logscale = 0_pInt ! linear/logaritmic time inc flag
logical :: followFormerTrajectory = . true . , & ! follow trajectory of former loadcase
velGradApplied = . false . ! decide wether velocity gradient or fdot is given
logical , dimension ( 3 , 3 ) :: maskDeformation = . false . , & ! mask of deformation boundary conditions
maskStress = . false . ! mask of stress boundary conditions
logical , dimension ( 9 ) :: maskStressVector = . false . ! linear mask of boundary conditions
end type
type ( bc_type ) , allocatable , dimension ( : ) :: bc
!--------------------------------------------------------------------------------------------------
! variables storing information from geom file
real ( pReal ) :: wgt
real ( pReal ) , dimension ( 3 ) :: geomdim = 0.0_pReal ! physical dimension of volume element per direction
integer ( pInt ) :: Npoints , & ! number of Fourier points
homog ! homogenization scheme used
integer ( pInt ) , dimension ( 3 ) :: res = 1_pInt ! resolution (number of Fourier points) in each direction
integer ( pInt ) :: res1_red ! to store res(1)/2 +1
!--------------------------------------------------------------------------------------------------
! stress, stiffness and compliance average etc.
2012-06-05 22:04:20 +05:30
real ( pReal ) , dimension ( 3 , 3 ) :: P_av = 0.0_pReal , &
2012-03-30 01:24:31 +05:30
F_aim = math_I3 , F_aim_lastInc = math_I3 , lambda_av , &
2012-03-07 23:07:40 +05:30
mask_stress , mask_defgrad , deltaF , F_star_av , &
F_aim_lab ! quantities rotated to other coordinate system
2012-06-05 22:04:20 +05:30
real ( pReal ) , dimension ( 3 , 3 , 3 , 3 ) :: C_inc0 , C = 0.0_pReal , S_inc0 , S_lastInc , C_lastInc ! stiffness and compliance
2012-03-07 23:07:40 +05:30
real ( pReal ) , dimension ( 6 ) :: sigma ! cauchy stress
real ( pReal ) , dimension ( 6 , 6 ) :: dsde ! small strain stiffness
real ( pReal ) , dimension ( 9 , 9 ) :: s_prev99 , c_prev99 ! compliance and stiffness in matrix notation
real ( pReal ) , dimension ( : , : ) , allocatable :: s_reduced , c_reduced ! reduced compliance and stiffness (only for stress BC)
integer ( pInt ) :: size_reduced = 0_pInt ! number of stress BCs
!--------------------------------------------------------------------------------------------------
! pointwise data
type ( C_PTR ) :: tensorField ! fields in real an fourier space
2012-06-05 22:04:20 +05:30
real ( pReal ) , dimension ( : , : , : , : , : ) , pointer :: tau_real , F_real ! fields in real space (pointer)
complex ( pReal ) , dimension ( : , : , : , : , : ) , pointer :: tau_fourier , F_fourier ! fields in fourier space (pointer)
real ( pReal ) , dimension ( : , : , : , : , : ) , allocatable :: F_lastInc , lambda , P , F_star
2012-03-31 18:11:46 +05:30
real ( pReal ) , dimension ( : , : , : , : , : , : , : ) , allocatable :: dPdF
2012-03-07 23:07:40 +05:30
real ( pReal ) , dimension ( : , : , : , : ) , allocatable :: coordinates
real ( pReal ) , dimension ( : , : , : ) , allocatable :: temperature
!--------------------------------------------------------------------------------------------------
! variables storing information for spectral method and FFTW
2012-06-05 22:04:20 +05:30
type ( C_PTR ) :: plan_correction , plan_tau ! plans for fftw
2012-03-07 23:07:40 +05:30
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
integer ( pInt ) , dimension ( 3 ) :: k_s
!--------------------------------------------------------------------------------------------------
! loop variables, convergence etc.
real ( pReal ) :: time = 0.0_pReal , time0 = 0.0_pReal , timeinc = 1.0_pReal , timeinc_old = 0.0_pReal ! elapsed time, begin of interval, time interval
2012-03-31 18:11:46 +05:30
real ( pReal ) :: guessmode , err_stress , err_stress_tol , err_f , err_p , err_crit , &
2012-06-05 22:04:20 +05:30
err_f_point , err_p_point , err_nl , err_nl_point
2012-03-07 23:07:40 +05:30
real ( pReal ) , dimension ( 3 , 3 ) , parameter :: ones = 1.0_pReal , zeroes = 0.0_pReal
complex ( pReal ) , dimension ( 3 , 3 ) :: temp33_Complex
real ( pReal ) , dimension ( 3 , 3 ) :: temp33_Real
integer ( pInt ) :: i , j , k , l , u , v , w , errorID = 0_pInt , ierr
2012-03-31 18:11:46 +05:30
integer ( pInt ) :: N_Loadcases , loadcase , inc , iter , ielem , CPFEM_mode , guesses , guessmax = 10_pInt , &
2012-03-07 23:07:40 +05:30
totalIncsCounter = 0_pInt , notConvergedCounter = 0_pInt , convergedCounter = 0_pInt
2012-03-31 18:11:46 +05:30
logical :: errmatinv , callCPFEM
2012-03-07 23:07:40 +05:30
character ( len = 6 ) :: loadcase_string
!--------------------------------------------------------------------------------------------------
!variables controlling debugging
logical :: debugGeneral , debugDivergence , debugRestart , debugFFTW
!##################################################################################################
! reading of information from load case file and geometry file
!##################################################################################################
open ( 6 , encoding = 'UTF-8' )
call DAMASK_interface_init
2012-04-11 22:58:08 +05:30
write ( 6 , '(a)' ) ''
write ( 6 , '(a)' ) ' <<<+- DAMASK_spectral_AL init -+>>>'
write ( 6 , '(a)' ) ' $Id$'
2012-03-07 23:07:40 +05:30
#include "compilation_info.f90"
2012-04-11 22:58:08 +05:30
write ( 6 , '(a)' ) ' Working Directory: ' , trim ( getSolverWorkingDirectoryName ( ) )
write ( 6 , '(a)' ) ' Solver Job Name: ' , trim ( getSolverJobName ( ) )
write ( 6 , '(a)' ) ''
2012-03-07 23:07:40 +05:30
!--------------------------------------------------------------------------------------------------
! reading the load case file and allocate data structure containing load cases
2012-04-11 22:58:08 +05:30
call IO_open_file ( myUnit , getLoadcaseName ( ) )
2012-03-07 23:07:40 +05:30
rewind ( myUnit )
do
read ( myUnit , '(a1024)' , END = 100 ) line
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
positions = IO_stringPos ( line , maxNchunksLoadcase )
do i = 1_pInt , maxNchunksLoadcase , 1_pInt ! reading compulsory parameters for loadcase
select case ( IO_lc ( IO_stringValue ( line , positions , i ) ) )
case ( 'l' , 'velocitygrad' , 'velgrad' , 'velocitygradient' )
N_l = N_l + 1_pInt
case ( 'fdot' , 'dotf' )
N_Fdot = N_Fdot + 1_pInt
case ( 't' , 'time' , 'delta' )
N_t = N_t + 1_pInt
case ( 'n' , 'incs' , 'increments' , 'steps' , 'logincs' , 'logsteps' )
N_n = N_n + 1_pInt
end select
enddo ! count all identifiers to allocate memory and do sanity check
enddo
100 N_Loadcases = N_n
if ( ( N_l + N_Fdot / = N_n ) . or . ( N_n / = N_t ) ) & ! sanity check
2012-04-11 22:58:08 +05:30
call IO_error ( error_ID = 837_pInt , ext_msg = trim ( getLoadcaseName ( ) ) ) ! error message for incomplete loadcase
2012-03-07 23:07:40 +05:30
allocate ( bc ( N_Loadcases ) )
!--------------------------------------------------------------------------------------------------
! reading the load case and assign values to the allocated data structure
rewind ( myUnit )
loadcase = 0_pInt
do
read ( myUnit , '(a1024)' , END = 101 ) line
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
loadcase = loadcase + 1_pInt
positions = IO_stringPos ( line , maxNchunksLoadcase )
do j = 1_pInt , maxNchunksLoadcase
select case ( IO_lc ( IO_stringValue ( line , positions , j ) ) )
case ( 'fdot' , 'dotf' , 'l' , 'velocitygrad' , 'velgrad' , 'velocitygradient' ) ! assign values for the deformation BC matrix
bc ( loadcase ) % velGradApplied = &
( IO_lc ( IO_stringValue ( line , positions , j ) ) == 'l' . or . & ! in case of given L, set flag to true
IO_lc ( IO_stringValue ( line , positions , j ) ) == 'velocitygrad' . or . &
IO_lc ( IO_stringValue ( line , positions , j ) ) == 'velgrad' . or . &
IO_lc ( IO_stringValue ( line , positions , j ) ) == 'velocitygradient' )
temp_valueVector = 0.0_pReal
temp_maskVector = . false .
forall ( k = 1_pInt : 9_pInt ) temp_maskVector ( k ) = IO_stringValue ( line , positions , j + k ) / = '*'
do k = 1_pInt , 9_pInt
if ( temp_maskVector ( k ) ) temp_valueVector ( k ) = IO_floatValue ( line , positions , j + k )
enddo
bc ( loadcase ) % maskDeformation = transpose ( reshape ( temp_maskVector , [ 3 , 3 ] ) )
bc ( loadcase ) % deformation = math_plain9to33 ( temp_valueVector )
case ( 'p' , 'pk1' , 'piolakirchhoff' , 'stress' )
temp_valueVector = 0.0_pReal
forall ( k = 1_pInt : 9_pInt ) bc ( loadcase ) % maskStressVector ( k ) = &
IO_stringValue ( line , positions , j + k ) / = '*'
do k = 1_pInt , 9_pInt
if ( bc ( loadcase ) % maskStressVector ( k ) ) temp_valueVector ( k ) = &
2012-03-20 23:31:31 +05:30
IO_floatValue ( line , positions , j + k ) ! assign values for the bc(loadcase)%P matrix
2012-03-07 23:07:40 +05:30
enddo
bc ( loadcase ) % maskStress = transpose ( reshape ( bc ( loadcase ) % maskStressVector , [ 3 , 3 ] ) )
2012-03-20 23:31:31 +05:30
bc ( loadcase ) % P = math_plain9to33 ( temp_valueVector )
2012-03-07 23:07:40 +05:30
case ( 't' , 'time' , 'delta' ) ! increment time
bc ( loadcase ) % time = IO_floatValue ( line , positions , j + 1_pInt )
case ( 'temp' , 'temperature' ) ! starting temperature
bc ( loadcase ) % temperature = IO_floatValue ( line , positions , j + 1_pInt )
case ( 'n' , 'incs' , 'increments' , 'steps' ) ! number of increments
bc ( loadcase ) % incs = IO_intValue ( line , positions , j + 1_pInt )
case ( 'logincs' , 'logincrements' , 'logsteps' ) ! number of increments (switch to log time scaling)
bc ( loadcase ) % incs = IO_intValue ( line , positions , j + 1_pInt )
bc ( loadcase ) % logscale = 1_pInt
case ( 'f' , 'freq' , 'frequency' , 'outputfreq' ) ! frequency of result writings
bc ( loadcase ) % outputfrequency = IO_intValue ( line , positions , j + 1_pInt )
case ( 'r' , 'restart' , 'restartwrite' ) ! frequency of writing restart information
bc ( loadcase ) % restartfrequency = max ( 0_pInt , IO_intValue ( line , positions , j + 1_pInt ) )
case ( 'guessreset' , 'dropguessing' )
bc ( loadcase ) % followFormerTrajectory = . false . ! do not continue to predict deformation along former trajectory
case ( 'euler' ) ! rotation of loadcase given in euler angles
u = 0_pInt ! assuming values given in radians
v = 1_pInt ! assuming keyword indicating degree/radians
select case ( IO_lc ( IO_stringValue ( line , positions , j + 1_pInt ) ) )
case ( 'deg' , 'degree' )
u = 1_pInt ! for conversion from degree to radian
case ( 'rad' , 'radian' )
case default
v = 0_pInt ! immediately reading in angles, assuming radians
end select
forall ( k = 1_pInt : 3_pInt ) temp33_Real ( k , 1 ) = &
IO_floatValue ( line , positions , j + v + k ) * real ( u , pReal ) * inRad
bc ( loadcase ) % rotation = math_EulerToR ( temp33_Real ( : , 1 ) )
case ( 'rotation' , 'rot' ) ! assign values for the rotation of loadcase matrix
temp_valueVector = 0.0_pReal
forall ( k = 1_pInt : 9_pInt ) temp_valueVector ( k ) = IO_floatValue ( line , positions , j + k )
bc ( loadcase ) % rotation = math_plain9to33 ( temp_valueVector )
end select
enddo ; enddo
101 close ( myUnit )
!-------------------------------------------------------------------------------------------------- ToDo: if temperature at CPFEM is treated properly, move this up immediately after interface init
! initialization of all related DAMASK modules (e.g. mesh.f90 reads in geometry)
call CPFEM_initAll ( bc ( 1 ) % temperature , 1_pInt , 1_pInt )
!--------------------------------------------------------------------------------------------------
2012-04-11 22:58:08 +05:30
! get resolution, dimension, homogenization and variables derived from resolution
2012-04-12 13:33:08 +05:30
res = mesh_spectral_getResolution ( )
geomdim = mesh_spectral_getDimension ( )
homog = mesh_spectral_getHomogenization ( )
2012-06-05 22:04:20 +05:30
2012-03-07 23:07:40 +05:30
res1_red = res ( 1 ) / 2_pInt + 1_pInt ! size of complex array in first dimension (c2r, r2c)
Npoints = res ( 1 ) * res ( 2 ) * res ( 3 )
wgt = 1.0_pReal / real ( Npoints , pReal )
2012-04-11 22:58:08 +05:30
2012-03-07 23:07:40 +05:30
!--------------------------------------------------------------------------------------------------
! output of geometry
2012-04-11 22:58:08 +05:30
write ( 6 , '(a)' ) ''
write ( 6 , '(a)' ) '#############################################################'
write ( 6 , '(a)' ) 'DAMASK spectral_AL:'
write ( 6 , '(a)' ) 'The AL spectral method boundary value problem solver for'
write ( 6 , '(a)' ) 'the Duesseldorf Advanced Material Simulation Kit'
write ( 6 , '(a)' ) '#############################################################'
write ( 6 , '(a)' ) 'geometry file: ' , trim ( getModelName ( ) ) / / InputFileExtension
write ( 6 , '(a)' ) '============================================================='
write ( 6 , '(a,3(i12 ))' ) 'resolution a b c:' , res
write ( 6 , '(a,3(f12.5))' ) 'dimension x y z:' , geomdim
write ( 6 , '(a,i5)' ) 'homogenization: ' , homog
write ( 6 , '(a)' ) '#############################################################'
write ( 6 , '(a)' ) 'loadcase file: ' , trim ( getLoadcaseName ( ) )
2012-03-07 23:07:40 +05:30
!--------------------------------------------------------------------------------------------------
! consistency checks and output of load case
bc ( 1 ) % followFormerTrajectory = . false . ! cannot guess along trajectory for first inc of first loadcase
do loadcase = 1_pInt , N_Loadcases
write ( loadcase_string , '(i6)' ) loadcase
2012-04-11 22:58:08 +05:30
write ( 6 , '(a)' ) '============================================================='
write ( 6 , '(a,i6)' ) 'loadcase: ' , loadcase
2012-03-07 23:07:40 +05:30
2012-04-11 22:58:08 +05:30
if ( . not . bc ( loadcase ) % followFormerTrajectory ) write ( 6 , '(a)' ) 'drop guessing along trajectory'
2012-03-07 23:07:40 +05:30
if ( bc ( loadcase ) % velGradApplied ) then
do j = 1_pInt , 3_pInt
if ( any ( bc ( loadcase ) % maskDeformation ( j , 1 : 3 ) . eqv . . true . ) . and . &
any ( bc ( loadcase ) % maskDeformation ( j , 1 : 3 ) . eqv . . false . ) ) errorID = 832_pInt ! each row should be either fully or not at all defined
enddo
2012-04-11 22:58:08 +05:30
write ( 6 , '(a)' ) 'velocity gradient:'
2012-03-07 23:07:40 +05:30
else
2012-04-11 22:58:08 +05:30
write ( 6 , '(a)' ) 'deformation gradient rate:'
2012-03-07 23:07:40 +05:30
endif
2012-04-11 22:58:08 +05:30
write ( 6 , '(3(3(f12.7,1x)/))' , advance = 'no' ) merge ( math_transpose33 ( bc ( loadcase ) % deformation ) , &
2012-03-07 23:07:40 +05:30
reshape ( spread ( DAMASK_NaN , 1 , 9 ) , [ 3 , 3 ] ) , transpose ( bc ( loadcase ) % maskDeformation ) )
2012-04-11 22:58:08 +05:30
write ( 6 , '(a,/,3(3(f12.7,1x)/))' , advance = 'no' ) 'stress / GPa:' , &
2012-03-20 23:31:31 +05:30
1e-9_pReal * merge ( math_transpose33 ( bc ( loadcase ) % P ) , &
2012-03-07 23:07:40 +05:30
reshape ( spread ( DAMASK_NaN , 1 , 9 ) , [ 3 , 3 ] ) , transpose ( bc ( loadcase ) % maskStress ) )
if ( any ( bc ( loadcase ) % rotation / = math_I3 ) ) &
write ( * , '(a,/,3(3(f12.7,1x)/))' , advance = 'no' ) ' rotation of loadframe:' , &
math_transpose33 ( bc ( loadcase ) % rotation )
2012-04-11 22:58:08 +05:30
write ( 6 , '(a,f12.6)' ) 'temperature:' , bc ( loadcase ) % temperature
write ( 6 , '(a,f12.6)' ) 'time: ' , bc ( loadcase ) % time
write ( 6 , '(a,i5)' ) 'increments: ' , bc ( loadcase ) % incs
write ( 6 , '(a,i5)' ) 'output frequency: ' , bc ( loadcase ) % outputfrequency
write ( 6 , '(a,i5)' ) 'restart frequency: ' , bc ( loadcase ) % restartfrequency
2012-03-07 23:07:40 +05:30
if ( any ( bc ( loadcase ) % maskStress . eqv . bc ( loadcase ) % maskDeformation ) ) errorID = 831_pInt ! exclusive or masking only
if ( any ( bc ( loadcase ) % maskStress . and . transpose ( bc ( loadcase ) % maskStress ) . and . &
reshape ( [ . false . , . true . , . true . , . true . , . false . , . true . , . true . , . true . , . false . ] , [ 3 , 3 ] ) ) ) &
errorID = 838_pInt ! no rotation is allowed by stress BC
if ( any ( abs ( math_mul33x33 ( bc ( loadcase ) % rotation , math_transpose33 ( bc ( loadcase ) % rotation ) ) &
- math_I3 ) > reshape ( spread ( rotation_tol , 1 , 9 ) , [ 3 , 3 ] ) ) &
. or . abs ( math_det33 ( bc ( loadcase ) % rotation ) ) > 1.0_pReal + rotation_tol ) &
errorID = 846_pInt ! given rotation matrix contains strain
if ( bc ( loadcase ) % time < 0.0_pReal ) errorID = 834_pInt ! negative time increment
if ( bc ( loadcase ) % incs < 1_pInt ) errorID = 835_pInt ! non-positive incs count
if ( bc ( loadcase ) % outputfrequency < 1_pInt ) errorID = 836_pInt ! non-positive result frequency
if ( errorID > 0_pInt ) call IO_error ( error_ID = errorID , ext_msg = loadcase_string )
enddo
!--------------------------------------------------------------------------------------------------
! debugging parameters
debugRestart = iand ( debug_spectral , debug_spectralRestart ) > 0_pInt
debugFFTW = iand ( debug_spectral , debug_spectralFFTW ) > 0_pInt
debugGeneral = . true .
!##################################################################################################
! initialization
!##################################################################################################
!--------------------------------------------------------------------------------------------------
! allocate more memory
2012-03-31 18:11:46 +05:30
allocate ( P ( res ( 1 ) , res ( 2 ) , res ( 3 ) , 3 , 3 ) , source = 0.0_pReal )
allocate ( dPdF ( res ( 1 ) , res ( 2 ) , res ( 3 ) , 3 , 3 , 3 , 3 ) , source = 0.0_pReal )
2012-03-07 23:07:40 +05:30
allocate ( F_star ( 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 )
allocate ( lambda ( res ( 1 ) , res ( 2 ) , res ( 3 ) , 3 , 3 ) , source = 0.0_pReal )
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 )
allocate ( temperature ( res ( 1 ) , res ( 2 ) , res ( 3 ) ) , source = bc ( 1 ) % temperature ) ! start out isothermally
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-06-05 22:04:20 +05:30
call c_f_pointer ( tensorField , tau_real , [ res ( 1 ) + 2_pInt , res ( 2 ) , res ( 3 ) , 3 , 3 ] ) ! place a pointer for the real representation
2012-03-07 23:07:40 +05:30
call c_f_pointer ( tensorField , F_real , [ res ( 1 ) + 2_pInt , res ( 2 ) , res ( 3 ) , 3 , 3 ] ) ! place a pointer for the real representation
2012-06-05 22:04:20 +05:30
call c_f_pointer ( tensorField , tau_fourier , [ res1_red , res ( 2 ) , res ( 3 ) , 3 , 3 ] ) ! place a pointer for the complex representation
2012-03-07 23:07:40 +05:30
call c_f_pointer ( tensorField , F_fourier , [ res1_red , res ( 2 ) , res ( 3 ) , 3 , 3 ] ) ! place a pointer for the complex representation
!--------------------------------------------------------------------------------------------------
! 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
2012-04-11 22:58:08 +05:30
!$ 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
2012-03-07 23:07:40 +05:30
call fftw_set_timelimit ( fftw_timelimit ) ! set timelimit for plan creation
!--------------------------------------------------------------------------------------------------
! creating plans
2012-06-05 22:04:20 +05:30
plan_tau = fftw_plan_many_dft_r2c ( 3 , [ res ( 3 ) , res ( 2 ) , res ( 1 ) ] , 9 , & ! dimensions , length in each dimension in reversed order
tau_real , [ res ( 3 ) , res ( 2 ) , res ( 1 ) + 2_pInt ] , & ! input data , physical length in each dimension in reversed order
2012-03-07 23:07:40 +05:30
1 , res ( 3 ) * res ( 2 ) * ( res ( 1 ) + 2_pInt ) , & ! striding , product of physical lenght in the 3 dimensions
2012-06-05 22:04:20 +05:30
tau_fourier , [ res ( 3 ) , res ( 2 ) , res1_red ] , &
2012-03-07 23:07:40 +05:30
1 , res ( 3 ) * res ( 2 ) * res1_red , fftw_planner_flag )
plan_correction = fftw_plan_many_dft_c2r ( 3 , [ res ( 3 ) , res ( 2 ) , res ( 1 ) ] , 9 , &
F_fourier , [ res ( 3 ) , res ( 2 ) , res1_red ] , &
1 , res ( 3 ) * res ( 2 ) * res1_red , &
F_real , [ res ( 3 ) , res ( 2 ) , res ( 1 ) + 2_pInt ] , &
1 , res ( 3 ) * res ( 2 ) * ( res ( 1 ) + 2_pInt ) , fftw_planner_flag )
2012-04-11 22:58:08 +05:30
if ( debugGeneral ) write ( 6 , '(a)' ) 'FFTW initialized'
2012-05-29 20:38:18 +05:30
!--------------------------------------------------------------------------------------------------
! init fields to no deformation
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
2012-06-05 22:04:20 +05:30
F_star ( i , j , k , 1 : 3 , 1 : 3 ) = math_I3
F_lastInc ( i , j , k , 1 : 3 , 1 : 3 ) = math_I3
2012-05-29 20:38:18 +05:30
coordinates ( i , j , k , 1 : 3 ) = geomdim / real ( res * [ i , j , k ] , pReal ) - geomdim / real ( 2_pInt * res , pReal )
call CPFEM_general ( 3_pInt , coordinates ( i , j , k , 1 : 3 ) , math_I3 , math_I3 , temperature ( i , j , k ) , &
0.0_pReal , ielem , 1_pInt , sigma , dsde , P ( i , j , k , 1 : 3 , 1 : 3 ) , dPdF ( i , j , k , 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 ) )
enddo ; enddo ; enddo
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 ( 2_pInt , coordinates ( i , j , k , 1 : 3 ) , math_I3 , math_I3 , temperature ( i , j , k ) , &
0.0_pReal , ielem , 1_pInt , sigma , dsde , P ( i , j , k , 1 : 3 , 1 : 3 ) , dPdF ( i , j , k , 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 ) )
C = C + dPdF ( i , j , k , 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 )
enddo ; enddo ; enddo
2012-06-05 22:04:20 +05:30
C_inc0 = C * wgt * 3.0_pReal ! linear reference material stiffness
2012-03-07 23:07:40 +05:30
!--------------------------------------------------------------------------------------------------
! calculation of discrete angular frequencies, ordered as in FFTW (wrap around) and remove the given highest frequencies
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 ) / geomdim
enddo ; enddo ; enddo
!--------------------------------------------------------------------------------------------------
! 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 , u = 1_pInt : 3_pInt ) &
xiDyad ( l , u ) = xi ( l , i , j , k ) * xi ( u , i , j , k )
forall ( l = 1_pInt : 3_pInt , u = 1_pInt : 3_pInt ) &
temp33_Real ( l , u ) = sum ( C_inc0 ( l , 1 : 3 , u , 1 : 3 ) * xiDyad )
temp33_Real = math_inv33 ( temp33_Real )
forall ( l = 1_pInt : 3_pInt , u = 1_pInt : 3_pInt , v = 1_pInt : 3_pInt , w = 1_pInt : 3_pInt ) &
gamma_hat ( i , j , k , l , u , v , w ) = temp33_Real ( l , v ) * xiDyad ( u , w )
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
!--------------------------------------------------------------------------------------------------
! write header of output file
open ( 538 , file = trim ( getSolverWorkingDirectoryName ( ) ) / / trim ( getSolverJobName ( ) ) &
/ / '.spectralOut' , form = 'UNFORMATTED' , status = 'REPLACE' )
write ( 538 ) 'load' , trim ( getLoadcaseName ( ) )
write ( 538 ) 'workingdir' , trim ( getSolverWorkingDirectoryName ( ) )
write ( 538 ) 'geometry' , trim ( getSolverJobName ( ) ) / / InputFileExtension
write ( 538 ) 'resolution' , res
write ( 538 ) 'dimension' , geomdim
write ( 538 ) 'materialpoint_sizeResults' , materialpoint_sizeResults
write ( 538 ) 'loadcases' , N_Loadcases
write ( 538 ) 'frequencies' , bc ( 1 : N_Loadcases ) % outputfrequency ! one entry per loadcase
write ( 538 ) 'times' , bc ( 1 : N_Loadcases ) % time ! one entry per loadcase
write ( 538 ) 'logscales' , bc ( 1 : N_Loadcases ) % logscale
write ( 538 ) 'increments' , bc ( 1 : N_Loadcases ) % incs ! one entry per loadcase
write ( 538 ) 'startingIncrement' , restartInc - 1_pInt ! start with writing out the previous inc
write ( 538 ) 'eoh' ! end of header
write ( 538 ) materialpoint_results ( 1_pInt : materialpoint_sizeResults , 1 , 1_pInt : Npoints ) ! initial (non-deformed or read-in) results
2012-04-11 22:58:08 +05:30
if ( debugGeneral ) write ( 6 , '(a)' ) 'Header of result file written out'
2012-03-07 23:07:40 +05:30
!##################################################################################################
! Loop over loadcases defined in the loadcase file
!##################################################################################################
do loadcase = 1_pInt , N_Loadcases
time0 = time ! loadcase start time
if ( bc ( loadcase ) % followFormerTrajectory . and . &
( restartInc < totalIncsCounter . or . &
restartInc > totalIncsCounter + bc ( loadcase ) % incs ) ) then ! continue to guess along former trajectory where applicable
guessmode = 1.0_pReal
else
guessmode = 0.0_pReal ! change of load case, homogeneous guess for the first inc
endif
!--------------------------------------------------------------------------------------------------
! arrays for mixed boundary conditions
mask_defgrad = merge ( ones , zeroes , bc ( loadcase ) % maskDeformation )
mask_stress = merge ( ones , zeroes , bc ( loadcase ) % maskStress )
size_reduced = int ( count ( bc ( loadcase ) % maskStressVector ) , pInt )
allocate ( c_reduced ( size_reduced , size_reduced ) , source = 0.0_pReal )
allocate ( s_reduced ( size_reduced , size_reduced ) , source = 0.0_pReal )
!##################################################################################################
! loop oper incs defined in input file for current loadcase
!##################################################################################################
do inc = 1_pInt , bc ( loadcase ) % incs
totalIncsCounter = totalIncsCounter + 1_pInt
if ( totalIncsCounter > = restartInc ) then ! do calculations (otherwise just forwarding)
!--------------------------------------------------------------------------------------------------
! forwarding time
timeinc_old = timeinc
if ( bc ( loadcase ) % logscale == 0_pInt ) then ! linear scale
timeinc = bc ( loadcase ) % time / bc ( loadcase ) % incs ! only valid for given linear time scale. will be overwritten later in case loglinear scale is used
else
if ( loadcase == 1_pInt ) then ! 1st loadcase of logarithmic scale
if ( inc == 1_pInt ) then ! 1st inc of 1st loadcase of logarithmic scale
timeinc = bc ( 1 ) % time * ( 2.0_pReal ** real ( 1_pInt - bc ( 1 ) % incs , pReal ) ) ! assume 1st inc is equal to 2nd
else ! not-1st inc of 1st loadcase of logarithmic scale
timeinc = bc ( 1 ) % time * ( 2.0_pReal ** real ( inc - 1_pInt - bc ( 1 ) % incs , pReal ) )
endif
else ! not-1st loadcase of logarithmic scale
timeinc = time0 * ( ( 1.0_pReal + bc ( loadcase ) % time / time0 ) ** ( real ( inc , pReal ) / &
real ( bc ( loadcase ) % incs , pReal ) ) &
- ( 1.0_pReal + bc ( loadcase ) % time / time0 ) ** ( real ( ( inc - 1_pInt ) , pReal ) / &
real ( bc ( loadcase ) % incs , pReal ) ) )
endif
endif
time = time + timeinc
if ( bc ( loadcase ) % velGradApplied ) then ! calculate deltaF from given L and current F
deltaF = timeinc * mask_defgrad * math_mul33x33 ( bc ( loadcase ) % deformation , F_aim )
else ! deltaF = fDot *timeinc where applicable
deltaF = timeinc * mask_defgrad * bc ( loadcase ) % deformation
endif
!--------------------------------------------------------------------------------------------------
! coordinates at beginning of inc
2012-03-30 01:24:31 +05:30
!call deformed_fft(res,geomdim,1.0_pReal,F_real(1:res(1),1:res(2),1:res(3),1:3,1:3),coordinates)! calculate current coordinates
2012-03-07 23:07:40 +05:30
!--------------------------------------------------------------------------------------------------
! 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
F_aim_lastInc = temp33_Real
2012-03-30 01:24:31 +05:30
2012-03-07 23:07:40 +05:30
!--------------------------------------------------------------------------------------------------
! update local deformation gradient
deltaF = math_rotate_backward33 ( deltaF , bc ( loadcase ) % rotation )
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 ) ; do i = 1_pInt , res ( 1 )
2012-06-05 22:04:20 +05:30
temp33_Real = F_star ( i , j , k , 1 : 3 , 1 : 3 )
F_star ( i , j , k , 1 : 3 , 1 : 3 ) = F_star ( i , j , k , 1 : 3 , 1 : 3 ) & ! decide if guessing along former trajectory or apply homogeneous addon
+ guessmode * ( F_star ( i , j , k , 1 : 3 , 1 : 3 ) - F_lastInc ( i , j , k , 1 : 3 , 1 : 3 ) ) & ! guessing...
2012-03-07 23:07:40 +05:30
* timeinc / timeinc_old &
+ ( 1.0_pReal - guessmode ) * deltaF ! if not guessing, use prescribed average deformation where applicable
F_lastInc ( i , j , k , 1 : 3 , 1 : 3 ) = temp33_Real
enddo ; enddo ; enddo
2012-03-31 18:11:46 +05:30
!--------------------------------------------------------------------------------------------------
! Initialize / Update lambda to useful value
2012-06-05 22:04:20 +05:30
P_av = P_av + math_mul3333xx33 ( C * wgt , F_aim - F_aim_lastInc )
2012-03-31 18:11:46 +05:30
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 ) ; do i = 1_pInt , res ( 1 )
2012-06-05 22:04:20 +05:30
lambda ( i , j , k , 1 : 3 , 1 : 3 ) = lambda ( i , j , k , 1 : 3 , 1 : 3 ) + math_mul3333xx33 ( C * wgt , F_aim - F_aim_lastInc )
2012-03-31 18:11:46 +05:30
enddo ; enddo ; enddo
2012-03-07 23:07:40 +05:30
!--------------------------------------------------------------------------------------------------
!Initialize pointwise data for AL scheme: ToDo: good choice?
2012-06-05 22:04:20 +05:30
F_star_av = 0.0_pReal
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 ) ; do i = 1_pInt , res ( 1 )
F_star_av = F_star_av + F_star ( i , j , k , 1 : 3 , 1 : 3 )
enddo ; enddo ; enddo
F_star_av = F_star_av * wgt
write ( * , '(a,/,3(3(f12.7,1x)/))' , advance = 'no' ) 'F* =' , &
math_transpose33 ( F_star_av )
!
2012-03-07 23:07:40 +05:30
!--------------------------------------------------------------------------------------------------
! calculate reduced compliance
if ( size_reduced > 0_pInt ) then ! calculate compliance in case stress BC is applied
C_lastInc = math_rotate_forward3333 ( C * wgt , bc ( loadcase ) % rotation ) ! calculate stiffness from former inc
c_prev99 = math_Plain3333to99 ( C_lastInc )
k = 0_pInt ! build reduced stiffness
do v = 1_pInt , 9_pInt
if ( bc ( loadcase ) % maskStressVector ( v ) ) then
k = k + 1_pInt
j = 0_pInt
do u = 1_pInt , 9_pInt
if ( bc ( loadcase ) % maskStressVector ( u ) ) then
j = j + 1_pInt
c_reduced ( k , j ) = c_prev99 ( v , u )
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 )
s_prev99 = 0.0_pReal ! build full compliance
k = 0_pInt
do v = 1_pInt , 9_pInt
if ( bc ( loadcase ) % maskStressVector ( v ) ) then
k = k + 1_pInt
j = 0_pInt
do u = 1_pInt , 9_pInt
if ( bc ( loadcase ) % maskStressVector ( u ) ) then
j = j + 1_pInt
s_prev99 ( v , u ) = s_reduced ( k , j )
endif ; enddo ; endif ; enddo
S_lastInc = ( math_Plain99to3333 ( s_prev99 ) )
endif
2012-06-05 22:04:20 +05:30
S_inc0 = math_invSym3333 ( C * wgt )
2012-03-07 23:07:40 +05:30
!--------------------------------------------------------------------------------------------------
! report begin of new increment
2012-04-11 22:58:08 +05:30
write ( 6 , '(a)' ) '##################################################################'
write ( 6 , '(A,I5.5,A,es12.5)' ) 'Increment ' , totalIncsCounter , ' Time ' , time
2012-03-07 23:07:40 +05:30
guessmode = 1.0_pReal ! keep guessing along former trajectory during same loadcase
CPFEM_mode = 1_pInt ! winding forward
iter = 0_pInt
2012-03-31 18:11:46 +05:30
err_crit = huge ( 1.0_pReal ) ! go into loop
callCPFEM = . true .
2012-05-29 20:38:18 +05:30
guessmax = 3
2012-03-31 18:11:46 +05:30
guesses = 0
2012-03-07 23:07:40 +05:30
!##################################################################################################
! convergence loop (looping over iterations)
!##################################################################################################
2012-03-31 18:11:46 +05:30
do while ( ( iter < itmax . and . ( err_crit > 1.0_pReal ) ) . or . iter < itmin )
2012-03-07 23:07:40 +05:30
iter = iter + 1_pInt
!--------------------------------------------------------------------------------------------------
! report begin of new iteration
2012-04-11 22:58:08 +05:30
write ( 6 , '(a)' ) ''
write ( 6 , '(a)' ) '=================================================================='
write ( 6 , '(5(a,i6.6))' ) 'Loadcase ' , loadcase , ' Increment ' , inc , '/' , bc ( loadcase ) % incs , &
2012-03-07 23:07:40 +05:30
' @ Iteration ' , iter , '/' , itmax
!--------------------------------------------------------------------------------------------------
! stress BC handling
if ( size_reduced > 0_pInt ) then ! calculate stress BC if applied
2012-06-05 22:04:20 +05:30
err_stress = maxval ( abs ( mask_stress * ( P_av - bc ( loadcase ) % P ) ) ) ! maximum deviaton (tensor norm not applicable)
F_aim = F_aim + math_mul3333xx33 ( S_lastInc , bc ( loadcase ) % P - P_av )
err_stress_tol = maxval ( abs ( P_av ) ) * err_stress_tolrel ! don't use any tensor norm because the comparison should be coherent
2012-03-07 23:07:40 +05:30
else
err_stress_tol = + huge ( 1.0_pReal )
endif
F_aim_lab = math_rotate_backward33 ( F_aim , bc ( loadcase ) % rotation )
2012-03-30 01:24:31 +05:30
write ( * , '(a,/,3(3(f12.7,1x)/))' , advance = 'no' ) 'F aim =' , &
2012-03-07 23:07:40 +05:30
math_transpose33 ( F_aim )
!--------------------------------------------------------------------------------------------------
! doing Fourier transform
2012-04-11 22:58:08 +05:30
write ( 6 , '(a)' ) '... spectral method ...............................................'
2012-06-05 22:04:20 +05:30
tau_real = 0.0_pReal
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 ) ; do i = 1_pInt , res ( 1 )
tau_real ( i , j , k , 1 : 3 , 1 : 3 ) = lambda ( i , j , k , 1 : 3 , 1 : 3 ) - math_mul3333xx33 ( C_inc0 , F_star ( i , j , k , 1 : 3 , 1 : 3 ) )
enddo ; enddo ; enddo
call fftw_execute_dft_r2c ( plan_tau , tau_real , tau_fourier )
tau_fourier ( res1_red , 1 : res ( 2 ) , 1 : res ( 3 ) , 1 : 3 , 1 : 3 ) &
2012-03-07 23:07:40 +05:30
= cmplx ( 0.0_pReal , 0.0_pReal , pReal )
2012-06-05 22:04:20 +05:30
tau_fourier ( 1 : res1_red , res ( 2 ) / 2_pInt + 1_pInt , 1 : res ( 3 ) , 1 : 3 , 1 : 3 ) &
2012-03-07 23:07:40 +05:30
= cmplx ( 0.0_pReal , 0.0_pReal , pReal )
if ( res ( 3 ) > 1_pInt ) &
2012-06-05 22:04:20 +05:30
tau_fourier ( 1 : res1_red , 1 : res ( 2 ) , res ( 3 ) / 2_pInt + 1_pInt , 1 : 3 , 1 : 3 ) &
2012-03-07 23:07:40 +05:30
= cmplx ( 0.0_pReal , 0.0_pReal , pReal )
2012-03-30 01:24:31 +05:30
!--------------------------------------------------------------------------------------------------
2012-03-07 23:07:40 +05:30
! using gamma operator to update F
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 , u = 1_pInt : 3_pInt ) &
xiDyad ( l , u ) = xi ( l , i , j , k ) * xi ( u , i , j , k )
forall ( l = 1_pInt : 3_pInt , u = 1_pInt : 3_pInt ) &
temp33_Real ( l , u ) = sum ( C_inc0 ( l , 1 : 3 , u , 1 : 3 ) * xiDyad )
temp33_Real = math_inv33 ( temp33_Real )
forall ( l = 1_pInt : 3_pInt , u = 1_pInt : 3_pInt , v = 1_pInt : 3_pInt , w = 1_pInt : 3_pInt ) &
gamma_hat ( 1 , 1 , 1 , l , u , v , w ) = temp33_Real ( l , v ) * xiDyad ( u , w )
forall ( l = 1_pInt : 3_pInt , u = 1_pInt : 3_pInt ) &
temp33_Complex ( l , u ) = sum ( gamma_hat ( 1 , 1 , 1 , l , u , 1 : 3 , 1 : 3 ) * &
2012-06-05 22:04:20 +05:30
tau_fourier ( i , j , k , 1 : 3 , 1 : 3 ) )
2012-03-07 23:07:40 +05:30
F_fourier ( i , j , k , 1 : 3 , 1 : 3 ) = - temp33_Complex
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 ( u = 1_pInt : 3_pInt , v = 1_pInt : 3_pInt ) &
temp33_Complex ( u , v ) = sum ( gamma_hat ( i , j , k , u , v , 1 : 3 , 1 : 3 ) * &
2012-06-05 22:04:20 +05:30
tau_fourier ( i , j , k , 1 : 3 , 1 : 3 ) )
2012-03-07 23:07:40 +05:30
F_fourier ( i , j , k , 1 : 3 , 1 : 3 ) = - temp33_Complex
enddo ; enddo ; enddo
endif
2012-06-05 22:04:20 +05:30
F_fourier ( 1 , 1 , 1 , 1 : 3 , 1 : 3 ) = cmplx ( ( F_aim_lab ) * real ( Npoints , pReal ) , 0.0_pReal , pReal )
2012-03-30 01:24:31 +05:30
2012-03-07 23:07:40 +05:30
!--------------------------------------------------------------------------------------------------
! doing inverse Fourier transform
call fftw_execute_dft_c2r ( plan_correction , F_fourier , F_real ) ! back transform of fluct deformation gradient
2012-06-05 22:04:20 +05:30
F_real ( 1 : res ( 1 ) , 1 : res ( 2 ) , 1 : res ( 3 ) , 1 : 3 , 1 : 3 ) = F_real ( 1 : res ( 1 ) , 1 : res ( 2 ) , 1 : res ( 3 ) , 1 : 3 , 1 : 3 ) * wgt
2012-03-07 23:07:40 +05:30
2012-06-05 22:04:20 +05:30
F_star_av = 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 )
F_star_av = F_star_av + F_star ( i , j , k , 1 : 3 , 1 : 3 )
temp33_real = temp33_real + F_real ( i , j , k , 1 : 3 , 1 : 3 )
enddo ; enddo ; enddo
F_star_av = F_star_av * wgt
write ( * , '(a,/,3(3(f12.7,1x)/))' , advance = 'no' ) 'F* =' , &
math_transpose33 ( F_star_av )
temp33_real = temp33_real * wgt
write ( * , '(a,/,3(3(f12.7,1x)/))' , advance = 'no' ) 'F =' , &
math_transpose33 ( temp33_real )
!--------------------------------------------------------------------------------------------------
write ( 6 , '(a)' ) '... calling CPFEM to update P(F*) and F*.........................'
err_nl_point = huge ( 1.0_pReal )
u = 0_pInt
do while ( err_nl_point > 1.0e6_pReal . or . u < 2_pInt )
u = u + 1_pInt
2012-03-31 18:11:46 +05:30
ielem = 0_pInt
2012-06-05 22:04:20 +05:30
err_nl_point = 0.0_pReal
err_nl = 0.0_pReal
2012-03-31 18:11:46 +05:30
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_star ( i , j , k , 1 : 3 , 1 : 3 ) , temperature ( i , j , k ) , timeinc , ielem , 1_pInt , &
sigma , dsde , P ( i , j , k , 1 : 3 , 1 : 3 ) , dPdF ( i , j , k , 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 ) )
enddo ; enddo ; enddo
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 ( CPFEM_mode , &
coordinates ( i , j , k , 1 : 3 ) , F_lastInc ( i , j , k , 1 : 3 , 1 : 3 ) , &
F_star ( i , j , k , 1 : 3 , 1 : 3 ) , temperature ( i , j , k ) , timeinc , ielem , 1_pInt , &
sigma , dsde , P ( i , j , k , 1 : 3 , 1 : 3 ) , dPdF ( i , j , k , 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 ) )
CPFEM_mode = 2_pInt ! winding forward
2012-06-05 22:04:20 +05:30
temp33_Real = P ( i , j , k , 1 : 3 , 1 : 3 ) - lambda ( i , j , k , 1 : 3 , 1 : 3 ) &
+ math_mul3333xx33 ( C_inc0 , F_star ( i , j , k , 1 : 3 , 1 : 3 ) - F_real ( i , j , k , 1 : 3 , 1 : 3 ) )
2012-03-31 18:11:46 +05:30
2012-06-05 22:04:20 +05:30
err_nl_point = max ( err_nl_point , maxval ( abs ( temp33_real ) ) )
err_nl = max ( err_nl , sqrt ( math_mul33xx33 ( temp33_real , temp33_real ) ) )
2012-03-31 18:11:46 +05:30
2012-06-05 22:04:20 +05:30
! temp33_Real = lambda(i,j,k,1:3,1:3) - P(i,j,k,1:3,1:3) &
! + math_mul3333xx33(C_inc0,F_real(i,j,k,1:3,1:3)- F_star(i,j,k,1:3,1:3))
! F_star(i,j,k,1:3,1:3) = F_star(i,j,k,1:3,1:3) + math_mul3333xx33(math_invSym3333(&
! C_inc0 + dPdF(i,j,k,1:3,1:3,1:3,1:3)), temp33_Real)
F_star ( i , j , k , 1 : 3 , 1 : 3 ) = F_star ( i , j , k , 1 : 3 , 1 : 3 ) - 0.1_pReal * math_mul3333xx33 ( S_inc0 , temp33_Real )
2012-03-31 18:11:46 +05:30
enddo ; enddo ; enddo
2012-06-05 22:04:20 +05:30
if ( u > 1 ) print * , 'comp wise error of nl eq last it:' , err_nl_point / 1.0e6_pReal
if ( u > 1 ) print * , 'norm error: last it ' , err_nl / 1.0e6_pReal
2012-03-30 01:24:31 +05:30
2012-06-05 22:04:20 +05:30
enddo
2012-03-31 18:11:46 +05:30
2012-06-05 22:04:20 +05:30
write ( 6 , '(a)' ) '... update λ..........................'
err_f = 0.0_pReal
err_f_point = 0.0_pReal
err_p = 0.0_pReal
err_p_point = 0.0_pReal
F_star_av = 0.0_pReal
lambda_av = 0.0_pReal
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 ) ; do i = 1_pInt , res ( 1 )
lambda ( i , j , k , 1 : 3 , 1 : 3 ) = lambda ( i , j , k , 1 : 3 , 1 : 3 ) + math_mul3333xx33 ( C_inc0 , F_real ( i , j , k , 1 : 3 , 1 : 3 ) &
- F_star ( i , j , k , 1 : 3 , 1 : 3 ) )
F_star_av = F_star_av + F_star ( i , j , k , 1 : 3 , 1 : 3 )
lambda_av = lambda_av + lambda ( i , j , k , 1 : 3 , 1 : 3 )
temp33_real = F_star ( i , j , k , 1 : 3 , 1 : 3 ) - F_real ( i , j , k , 1 : 3 , 1 : 3 )
err_f_point = max ( err_f_point , maxval ( abs ( temp33_real ) ) )
err_f = max ( err_f , sqrt ( math_mul33xx33 ( temp33_real , temp33_real ) ) )
enddo ; enddo ; enddo
ielem = 0_pInt
P_av = 0.0_pReal
C = 0.0_pReal
2012-03-07 23:07:40 +05:30
do k = 1_pInt , res ( 3 ) ; do j = 1_pInt , res ( 2 ) ; do i = 1_pInt , res ( 1 )
2012-06-05 22:04:20 +05:30
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_real ( i , j , k , 1 : 3 , 1 : 3 ) , temperature ( i , j , k ) , timeinc , ielem , 1_pInt , &
sigma , dsde , P ( i , j , k , 1 : 3 , 1 : 3 ) , dPdF ( i , j , k , 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 ) )
enddo ; enddo ; enddo
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 ( 2_pInt , &
coordinates ( i , j , k , 1 : 3 ) , F_lastInc ( i , j , k , 1 : 3 , 1 : 3 ) , &
F_real ( i , j , k , 1 : 3 , 1 : 3 ) , temperature ( i , j , k ) , timeinc , ielem , 1_pInt , &
sigma , dsde , P ( i , j , k , 1 : 3 , 1 : 3 ) , dPdF ( i , j , k , 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 ) )
P_av = P_av + P ( i , j , k , 1 : 3 , 1 : 3 )
temp33_real = lambda ( i , j , k , 1 : 3 , 1 : 3 ) - P ( i , j , k , 1 : 3 , 1 : 3 )
2012-03-31 18:11:46 +05:30
err_p_point = max ( err_p_point , maxval ( abs ( temp33_real ) ) )
2012-03-30 01:24:31 +05:30
err_p = max ( err_p , sqrt ( math_mul33xx33 ( temp33_real , temp33_real ) ) )
2012-06-05 22:04:20 +05:30
C = C + dPdF ( i , j , k , 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 )
2012-03-07 23:07:40 +05:30
enddo ; enddo ; enddo
2012-03-20 23:31:31 +05:30
2012-03-07 23:07:40 +05:30
F_star_av = F_star_av * wgt
2012-03-30 01:24:31 +05:30
write ( * , '(a,/,3(3(f12.7,1x)/))' , advance = 'no' ) 'F* =' , &
2012-03-07 23:07:40 +05:30
math_transpose33 ( F_star_av )
2012-06-05 22:04:20 +05:30
P_av = P_av * wgt
write ( * , '(a,/,3(3(es14.7,1x)/))' , advance = 'no' ) 'P(F) / GPa =' , &
math_transpose33 ( P_av ) / 1.e6_pReal
2012-03-30 01:24:31 +05:30
lambda_av = lambda_av * wgt
write ( * , '(a,/,3(3(es14.7,1x)/))' , advance = 'no' ) 'λ / GPa =' , &
math_transpose33 ( lambda_av ) / 1.e6_pReal
2012-03-07 23:07:40 +05:30
err_f = err_f / sqrt ( math_mul33xx33 ( F_star_av , F_star_av ) )
2012-06-05 22:04:20 +05:30
err_p = err_p / sqrt ( math_mul33xx33 ( P_av , P_av ) )
2012-03-31 18:11:46 +05:30
write ( 6 , '(a,es14.7,es14.7)' ) 'error F' , err_f / 1e-4 , err_f
write ( 6 , '(a,es14.7,es14.7)' ) 'error P' , err_p / 1e-3 , err_p
write ( 6 , '(a,es14.7,es14.7)' ) 'error stress = ' , err_stress / err_stress_tol , err_stress
write ( 6 , * ) ' '
write ( 6 , '(a,es14.7)' ) 'max abs err F' , err_f_point
write ( 6 , '(a,es14.7)' ) 'max abs err P' , err_p_point
2012-06-05 22:04:20 +05:30
err_crit = max ( err_p / 1e-3 , err_f / 1e-4 , err_stress / err_stress_tol )
2012-04-12 13:33:08 +05:30
write ( 6 , '(a,es14.7)' ) 'critical error' , err_crit
2012-03-31 18:11:46 +05:30
2012-03-07 23:07:40 +05:30
enddo ! end looping when convergency is achieved
2012-04-06 19:53:06 +05:30
write ( 6 , '(a)' ) ' '
write ( 6 , '(a)' ) '=================================================================='
2012-03-31 18:11:46 +05:30
if ( err_crit > 1.0_pReal ) then
2012-04-06 19:53:06 +05:30
write ( 6 , '(A,I5.5,A)' ) 'increment ' , totalIncsCounter , ' NOT converged'
2012-03-07 23:07:40 +05:30
notConvergedCounter = notConvergedCounter + 1_pInt
else
convergedCounter = convergedCounter + 1_pInt
2012-04-06 19:53:06 +05:30
write ( 6 , '(A,I5.5,A)' ) 'increment ' , totalIncsCounter , ' converged'
2012-03-07 23:07:40 +05:30
endif
if ( mod ( totalIncsCounter - 1_pInt , bc ( loadcase ) % outputfrequency ) == 0_pInt ) then ! at output frequency
2012-04-06 19:53:06 +05:30
write ( 6 , '(a)' ) ' '
write ( 6 , '(a)' ) '... writing results to file ......................................'
2012-03-07 23:07:40 +05:30
write ( 538 ) materialpoint_results ( 1_pInt : materialpoint_sizeResults , 1 , 1_pInt : Npoints ) ! write result to file
endif
endif ! end calculation/forwarding
enddo ! end looping over incs in current loadcase
deallocate ( c_reduced )
deallocate ( s_reduced )
enddo ! end looping over loadcases
2012-04-11 22:58:08 +05:30
write ( 6 , '(a)' ) ''
write ( 6 , '(a)' ) '##################################################################'
write ( 6 , '(i6.6,a,i6.6,a)' ) notConvergedCounter , ' out of ' , &
notConvergedCounter + convergedCounter , ' increments did not converge!'
2012-03-07 23:07:40 +05:30
close ( 538 )
2012-06-05 22:04:20 +05:30
call fftw_destroy_plan ( plan_tau ) ; call fftw_destroy_plan ( plan_correction )
2012-03-31 18:11:46 +05:30
call quit ( 1_pInt )
2012-04-12 13:33:08 +05:30
end program DAMASK_spectral_AL