2010-06-08 15:40:57 +05:30
!* $Id$
2010-06-08 15:38:15 +05:30
!********************************************************************
! Material subroutine for BVP solution using spectral method
!
! written by P. Eisenlohr,
! F. Roters,
! L. Hantcherli,
! W.A. Counts
! D.D. Tjahjanto
! C. Kords
! M. Diehl
! R. Lebensohn
!
! MPI fuer Eisenforschung, Duesseldorf
!
!********************************************************************
! Usage:
! - start program with mpie_spectral PathToMeshFile/NameOfMesh.mesh
! PathToLoadFile/NameOfLoadFile.load
! - PathToLoadFile will be the working directory
! - make sure the file "material.config" exists in the working
! directory
!********************************************************************
!
include "prec.f90" ! uses nothing else
MODULE mpie_interface
2010-06-10 20:21:10 +05:30
use prec , only : pInt , pReal
character ( len = 64 ) , parameter :: FEsolver = 'Spectral'
character ( len = 5 ) , parameter :: InputFileExtension = '.mesh'
2010-06-08 15:38:15 +05:30
CONTAINS
2010-06-10 20:21:10 +05:30
!********************************************************************
! initialize interface module
!
!********************************************************************
subroutine mpie_interface_init ( )
2010-06-08 15:38:15 +05:30
write ( 6 , * )
2010-06-08 15:40:57 +05:30
write ( 6 , * ) '<<<+- mpie_spectral init -+>>>'
write ( 6 , * ) '$Id$'
2010-06-08 15:38:15 +05:30
write ( 6 , * )
2010-06-10 20:21:10 +05:30
2010-06-08 15:38:15 +05:30
return
2010-07-02 19:40:36 +05:30
endsubroutine
2010-06-08 15:38:15 +05:30
2010-06-10 20:21:10 +05:30
!********************************************************************
! extract working directory from loadcase file
! possibly based on current working dir
!********************************************************************
2010-06-08 15:38:15 +05:30
function getSolverWorkingDirectoryName ( )
2010-06-10 20:21:10 +05:30
2010-06-08 15:38:15 +05:30
implicit none
2010-06-10 20:21:10 +05:30
character ( len = 1024 ) cwd , outname , getSolverWorkingDirectoryName
character ( len = * ) , parameter :: pathSep = achar ( 47 ) / / achar ( 92 ) ! /, \
2010-06-08 15:38:15 +05:30
call getarg ( 2 , outname ) ! path to loadFile
if ( scan ( outname , pathSep ) == 1 ) then ! absolute path given as command line argument
getSolverWorkingDirectoryName = outname ( 1 : scan ( outname , pathSep , back = . true . ) )
else
call getcwd ( cwd )
getSolverWorkingDirectoryName = trim ( cwd ) / / '/' / / outname ( 1 : scan ( outname , pathSep , back = . true . ) )
endif
2010-06-10 20:21:10 +05:30
getSolverWorkingDirectoryName = rectifyPath ( getSolverWorkingDirectoryName )
2010-06-08 15:38:15 +05:30
return
2010-06-10 20:21:10 +05:30
2010-07-02 19:40:36 +05:30
endfunction
2010-06-08 15:38:15 +05:30
2010-06-10 20:21:10 +05:30
!********************************************************************
! basename of meshfile from command line arguments
!
!********************************************************************
2010-06-08 15:38:15 +05:30
function getSolverJobName ( )
2010-06-10 20:21:10 +05:30
use prec , only : pInt
2010-06-08 15:38:15 +05:30
implicit none
2010-06-10 20:21:10 +05:30
character ( 1024 ) getSolverJobName , outName , cwd
2010-06-08 15:38:15 +05:30
character ( len = * ) , parameter :: pathSep = achar ( 47 ) / / achar ( 92 ) ! /, \
2010-06-10 20:21:10 +05:30
integer ( pInt ) posExt , posSep
getSolverJobName = ''
2010-06-08 15:38:15 +05:30
call getarg ( 1 , outName )
2010-06-10 20:21:10 +05:30
posExt = scan ( outName , '.' , back = . true . )
posSep = scan ( outName , pathSep , back = . true . )
2010-07-01 20:50:06 +05:30
if ( posExt < = posSep ) posExt = len_trim ( outName ) + 1 ! no extension present
getSolverJobName = outName ( 1 : posExt - 1 ) ! path to mesh file (excl. extension)
if ( scan ( getSolverJobName , pathSep ) / = 1 ) then ! relative path given as command line argument
2010-06-10 20:21:10 +05:30
call getcwd ( cwd )
2010-07-01 20:50:06 +05:30
getSolverJobName = rectifyPath ( trim ( cwd ) / / '/' / / getSolverJobName )
else
getSolverJobName = rectifyPath ( getSolverJobName )
2010-06-10 20:21:10 +05:30
endif
2010-07-01 20:50:06 +05:30
getSolverJobName = makeRelativePath ( getSolverWorkingDirectoryName ( ) , &
getSolverJobName )
2010-06-10 20:21:10 +05:30
return
2010-07-02 19:40:36 +05:30
endfunction
2010-06-10 20:21:10 +05:30
!********************************************************************
2010-07-01 20:50:06 +05:30
! relative path of loadcase from command line arguments
2010-06-10 20:21:10 +05:30
!
!********************************************************************
function getLoadcaseName ( )
use prec , only : pInt
implicit none
character ( len = 1024 ) getLoadcaseName , outName , cwd
character ( len = * ) , parameter :: pathSep = achar ( 47 ) / / achar ( 92 ) ! /, \
integer ( pInt ) posExt , posSep
2010-07-02 19:40:36 +05:30
posExt = 0 !not sure if its needed
2010-06-10 20:21:10 +05:30
call getarg ( 2 , getLoadcaseName )
2010-07-01 20:50:06 +05:30
posExt = scan ( getLoadcaseName , '.' , back = . true . )
posSep = scan ( getLoadcaseName , pathSep , back = . true . )
if ( posExt < = posSep ) getLoadcaseName = trim ( getLoadcaseName ) / / ( '.load' ) ! no extension present
2010-06-10 20:21:10 +05:30
if ( scan ( getLoadcaseName , pathSep ) / = 1 ) then ! relative path given as command line argument
call getcwd ( cwd )
getLoadcaseName = rectifyPath ( trim ( cwd ) / / '/' / / getLoadcaseName )
2010-07-01 20:50:06 +05:30
else
getLoadcaseName = rectifyPath ( getLoadcaseName )
2010-06-10 20:21:10 +05:30
endif
getLoadcaseName = makeRelativePath ( getSolverWorkingDirectoryName ( ) , &
getLoadcaseName )
return
2010-07-02 19:40:36 +05:30
endfunction
2010-06-08 15:38:15 +05:30
2010-06-10 14:20:04 +05:30
2010-06-10 20:21:10 +05:30
!********************************************************************
! remove ../ and ./ from path
!
!********************************************************************
function rectifyPath ( path )
use prec , only : pInt
implicit none
character ( len = * ) path
character ( len = len_trim ( path ) ) rectifyPath
integer ( pInt ) i , j , k , l
!remove ./ from path
l = len_trim ( path )
rectifyPath = path
2010-07-02 19:40:36 +05:30
do i = l , 2 , - 1
if ( rectifyPath ( i - 1 : i ) == './' . and . rectifyPath ( i - 2 : i - 2 ) / = '.' ) &
2010-06-10 14:20:04 +05:30
rectifyPath ( i - 1 : l ) = rectifyPath ( i + 1 : l ) / / ' '
2010-07-02 19:40:36 +05:30
enddo
2010-06-10 20:21:10 +05:30
!remove ../ and corresponding directory from rectifyPath
l = len_trim ( rectifyPath )
i = index ( rectifyPath ( i : l ) , '../' )
j = 0_pInt
do while ( i > j )
j = scan ( rectifyPath ( : i - 2 ) , '/' , back = . true . )
rectifyPath ( j + 1 : l ) = rectifyPath ( i + 3 : l ) / / repeat ( ' ' , 2 + i - j )
i = j + index ( rectifyPath ( j + 1 : l ) , '../' )
2010-07-02 19:40:36 +05:30
enddo
if ( len_trim ( rectifyPath ) == 0 ) rectifyPath = '/'
2010-06-10 20:21:10 +05:30
return
2010-07-02 19:40:36 +05:30
endfunction rectifyPath
2010-06-10 20:21:10 +05:30
!********************************************************************
! relative path from absolute a to absolute b
!
!********************************************************************
2010-06-10 14:20:04 +05:30
function makeRelativePath ( a , b )
2010-06-10 20:21:10 +05:30
use prec , only : pInt
2010-06-10 14:20:04 +05:30
implicit none
2010-06-10 20:21:10 +05:30
character ( len = * ) :: a , b
character ( len = 1024 ) :: makeRelativePath
integer ( pInt ) i , posLastCommonSlash , remainingSlashes
2010-06-10 14:20:04 +05:30
posLastCommonSlash = 0
remainingSlashes = 0
2010-06-10 20:21:10 +05:30
do i = 1 , min ( 1024 , len_trim ( a ) , len_trim ( b ) )
2010-06-10 14:20:04 +05:30
if ( a ( i : i ) / = b ( i : i ) ) exit
if ( a ( i : i ) == '/' ) posLastCommonSlash = i
enddo
do i = posLastCommonSlash + 1 , len_trim ( a )
if ( a ( i : i ) == '/' ) remainingSlashes = remainingSlashes + 1
enddo
2010-06-10 20:21:10 +05:30
makeRelativePath = repeat ( '../' , remainingSlashes ) / / b ( posLastCommonSlash + 1 : len_trim ( b ) )
return
2010-07-02 19:40:36 +05:30
endfunction makeRelativePath
2010-06-10 14:20:04 +05:30
2010-06-08 15:38:15 +05:30
END MODULE
2010-06-10 20:21:10 +05:30
include "IO.f90" ! uses prec
include "numerics.f90" ! uses prec, IO
include "math.f90" ! uses prec, numerics
include "debug.f90" ! uses prec, numerics
include "FEsolving.f90" ! uses prec, IO
include "mesh.f90" ! uses prec, math, IO, FEsolving
include "material.f90" ! uses prec, math, IO, mesh
include "lattice.f90" ! uses prec, math, IO, material
include "constitutive_phenopowerlaw.f90" ! uses prec, math, IO, lattice, material, debug
include "constitutive_j2.f90" ! uses prec, math, IO, lattice, material, debug
include "constitutive_dislotwin.f90" ! uses prec, math, IO, lattice, material, debug
include "constitutive_nonlocal.f90" ! uses prec, math, IO, lattice, material, debug
include "constitutive.f90" ! uses prec, IO, math, lattice, mesh, debug
include "crystallite.f90" ! uses prec, math, IO, numerics
include "homogenization_isostrain.f90" ! uses prec, math, IO,
include "homogenization_RGC.f90" ! uses prec, math, IO, numerics, mesh: added <<<updated 31.07.2009>>>
include "homogenization.f90" ! uses prec, math, IO, numerics
include "CPFEM.f90" ! uses prec, math, IO, numerics, debug, FEsolving, mesh, lattice, constitutive, crystallite
2010-06-08 15:38:15 +05:30
2010-06-10 20:21:10 +05:30
!********************************************************************
2010-06-08 15:38:15 +05:30
program mpie_spectral
2010-06-10 20:21:10 +05:30
!********************************************************************
2010-06-08 15:38:15 +05:30
use mpie_interface
2010-06-10 20:21:10 +05:30
use prec , only : pInt , pReal
use IO
2010-07-01 20:50:06 +05:30
! use math, only: math_I3,math_transpose3x3,math_Mandel66to3333
use math
2010-06-10 21:02:06 +05:30
use CPFEM , only : CPFEM_general
2010-06-10 20:21:10 +05:30
2010-06-08 15:38:15 +05:30
implicit none
2010-06-10 20:21:10 +05:30
2010-06-10 21:02:06 +05:30
real ( pReal ) , dimension ( : , : , : ) , allocatable :: bc_velocityGrad , &
bc_stress ! velocity gradient and stress BC
real ( pReal ) , dimension ( : ) , allocatable :: bc_timeIncrement ! length of increment
integer ( pInt ) , dimension ( : ) , allocatable :: bc_steps ! number of steps
logical , dimension ( : , : , : , : ) , allocatable :: bc_mask ! mask
real ( pReal ) temperature
real ( pReal ) , dimension ( 6 ) :: stress
real ( pReal ) , dimension ( 6 , 6 ) :: dsde
2010-06-08 15:38:15 +05:30
2010-06-10 20:21:10 +05:30
character ( len = 1024 ) path , line
2010-06-10 21:02:06 +05:30
logical , dimension ( 9 ) :: bc_maskvector
2010-06-25 17:01:05 +05:30
logical gotResolution , gotDimension , gotHomogenization
integer ( pInt ) , parameter :: maxNchunksInput = 24 ! 4 identifiers, 18 values for the matrices and 2 scalars
integer ( pInt ) , dimension ( 1 + maxNchunksInput * 2 ) :: posInput
integer ( pInt ) , parameter :: maxNchunksMesh = 7 ! 4 identifiers, 3 values
integer ( pInt ) , dimension ( 1 + 2 * maxNchunksMesh ) :: posMesh
2010-06-10 20:21:10 +05:30
real ( pReal ) , dimension ( 9 ) :: valuevector
2010-06-25 17:01:05 +05:30
integer ( pInt ) unit , N_l , N_s , N_t , N_n , N , i , j , k , l ! numbers of identifiers, loop variables
integer ( pInt ) a , b , c , e , homog
real ( pReal ) x , y , z
2010-07-02 19:40:36 +05:30
2010-07-01 20:50:06 +05:30
!-------------------------
!begin RL
!-------------------------
real ( pReal ) , dimension ( : ) , allocatable :: datafft
real ( pReal ) , dimension ( : , : , : , : , : ) , allocatable :: workfft , workfftim , sg , disgrad , defgradold
integer ( pInt ) , dimension ( 3 , 3 ) :: iudot , iscau
real ( pReal ) , dimension ( 3 , 3 ) :: disgradmacro , disgradmacroactual
real ( pReal ) , dimension ( 3 , 3 ) :: ddisgradmacro , ddisgradmacroacum , ddisgrad , ddisgradim
real ( pReal ) , dimension ( 3 , 3 ) :: defgrad0 , defgrad
2010-07-02 22:45:53 +05:30
real ( pReal ) , dimension ( 3 , 3 ) :: udot , scauchy , scauav , aux33 , xkdyad
2010-07-01 20:50:06 +05:30
2010-07-02 19:40:36 +05:30
integer ( pInt ) , dimension ( 3 ) :: nn
integer ( pInt ) , dimension ( 2 ) :: nn2
2010-07-05 17:03:48 +05:30
real ( pReal ) , dimension ( 3 ) :: delt , xk
2010-07-01 20:50:06 +05:30
real ( pReal ) , dimension ( 6 ) :: aux6
real ( pReal ) , dimension ( 3 , 3 , 3 , 3 ) :: c0 , s0 , g1
real ( pReal ) , dimension ( 6 , 6 ) :: c066 , s066
2010-07-02 19:40:36 +05:30
integer ( pInt ) itmax , jload , ielem , ii , jj , k1 , kxx , kyy , kzz , kx , ky , kz , idum , iter , imicro , m1 , n1 , p , q
2010-07-01 20:50:06 +05:30
real ( pReal ) prodnn , wgt , error , tdot , erre , errs , evm , svm , det , xknorm
logical errmatinv
!-------------------------
!end RL
!-------------------------
2010-06-25 17:01:05 +05:30
2010-06-10 20:21:10 +05:30
if ( IargC ( ) < 2 ) call IO_error ( 102 )
path = getLoadcaseName ( )
2010-06-10 21:02:06 +05:30
bc_maskvector = ''
2010-06-10 20:21:10 +05:30
unit = 234_pInt
N_l = 0_pInt
N_s = 0_pInt
N_t = 0_pInt
N_n = 0_pInt
2010-07-01 20:50:06 +05:30
print * , 'Loadcase: ' , trim ( path )
print * , 'Workingdir: ' , trim ( getSolverWorkingDirectoryName ( ) )
2010-06-10 20:21:10 +05:30
2010-07-02 19:40:36 +05:30
if ( . not . IO_open_file ( unit , path ) ) call IO_error ( 45 , ext_msg = path )
2010-06-10 20:21:10 +05:30
rewind ( unit )
do
2010-07-02 19:40:36 +05:30
read ( unit , '(a1024)' , END = 101 ) line
2010-06-10 20:21:10 +05:30
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
2010-06-25 17:01:05 +05:30
posInput = IO_stringPos ( line , maxNchunksInput )
do i = 1 , maxNchunksInput , 1
select case ( IO_lc ( IO_stringValue ( line , posInput , i ) ) )
2010-06-10 20:21:10 +05:30
case ( 'l' , 'velocitygrad' )
N_l = N_l + 1
case ( 's' , 'stress' )
N_s = N_s + 1
case ( 't' , 'time' , 'delta' )
N_t = N_t + 1
case ( 'n' , 'incs' , 'increments' , 'steps' )
N_n = N_n + 1
end select
enddo ! count all identifiers to allocate memory and do sanity check
if ( ( N_l / = N_s ) . or . ( N_s / = N_t ) . or . ( N_t / = N_n ) ) & ! sanity check
2010-07-02 19:40:36 +05:30
call IO_error ( 46 , ext_msg = path ) ! error message for incomplete input file
2010-06-10 20:21:10 +05:30
enddo
2010-06-10 14:20:04 +05:30
! allocate memory depending on lines in input file
2010-06-10 20:21:10 +05:30
101 N = N_l
2010-06-10 21:02:06 +05:30
allocate ( bc_velocityGrad ( 3 , 3 , N ) ) ; bc_velocityGrad = 0.0_pReal
allocate ( bc_stress ( 3 , 3 , N ) ) ; bc_stress = 0.0_pReal
allocate ( bc_mask ( 3 , 3 , 2 , N ) ) ; bc_mask = . false .
allocate ( bc_timeIncrement ( N ) ) ; bc_timeIncrement = 0.0_pReal
allocate ( bc_steps ( N ) ) ; bc_steps = 0_pInt
2010-06-10 20:21:10 +05:30
rewind ( unit )
j = 0_pInt
do
2010-07-02 19:40:36 +05:30
read ( unit , '(a1024)' , END = 200 ) line
2010-06-10 14:20:04 +05:30
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
2010-07-02 19:40:36 +05:30
j = j + 1
2010-06-25 17:01:05 +05:30
posInput = IO_stringPos ( line , maxNchunksInput )
do i = 1 , maxNchunksInput , 2
select case ( IO_lc ( IO_stringValue ( line , posInput , i ) ) )
2010-06-10 20:21:10 +05:30
case ( 'l' , 'velocitygrad' )
valuevector = 0.0_pReal
2010-06-25 17:01:05 +05:30
forall ( k = 1 : 9 ) bc_maskvector ( k ) = IO_stringValue ( line , posInput , i + k ) / = '#'
2010-06-10 20:21:10 +05:30
do k = 1 , 9
2010-06-25 17:01:05 +05:30
if ( bc_maskvector ( k ) ) valuevector ( k ) = IO_floatValue ( line , posInput , i + k ) ! assign values for the velocity gradient matrix
2010-06-10 20:21:10 +05:30
enddo
2010-06-10 21:02:06 +05:30
bc_mask ( : , : , 1 , j ) = reshape ( bc_maskvector , ( / 3 , 3 / ) )
bc_velocityGrad ( : , : , j ) = reshape ( valuevector , ( / 3 , 3 / ) )
2010-06-10 20:21:10 +05:30
case ( 's' , 'stress' )
valuevector = 0.0_pReal
2010-06-25 17:01:05 +05:30
forall ( k = 1 : 9 ) bc_maskvector ( k ) = IO_stringValue ( line , posInput , i + k ) / = '#'
2010-06-10 20:21:10 +05:30
do k = 1 , 9
2010-06-25 17:01:05 +05:30
if ( bc_maskvector ( k ) ) valuevector ( k ) = IO_floatValue ( line , posInput , i + k ) ! assign values for the bc_stress matrix
2010-06-10 20:21:10 +05:30
enddo
2010-06-10 21:02:06 +05:30
bc_mask ( : , : , 2 , j ) = reshape ( bc_maskvector , ( / 3 , 3 / ) )
bc_stress ( : , : , j ) = reshape ( valuevector , ( / 3 , 3 / ) )
2010-06-10 20:21:10 +05:30
case ( 't' , 'time' , 'delta' ) ! increment time
2010-06-25 17:01:05 +05:30
bc_timeIncrement ( j ) = IO_floatValue ( line , posInput , i + 1 )
2010-06-10 21:02:06 +05:30
case ( 'n' , 'incs' , 'increments' , 'steps' ) ! bc_steps
2010-06-25 17:01:05 +05:30
bc_steps ( j ) = IO_intValue ( line , posInput , i + 1 )
2010-06-10 20:21:10 +05:30
end select
2010-06-10 14:20:04 +05:30
enddo
2010-06-10 20:21:10 +05:30
enddo
200 close ( unit )
! consistency checks
do j = 1 , N
2010-06-10 21:02:06 +05:30
if ( any ( bc_mask ( : , : , 1 , j ) == bc_mask ( : , : , 2 , j ) ) ) &
call IO_error ( 47 , j ) ! bc_mask consistency
if ( any ( math_transpose3x3 ( bc_stress ( : , : , j ) ) + bc_stress ( : , : , j ) / = 2.0_pReal * bc_stress ( : , : , j ) ) ) &
call IO_error ( 48 , j ) ! bc_stress symmetry
2010-06-10 20:21:10 +05:30
2010-06-10 21:02:06 +05:30
print '(a,/,3(3(f12.6,x)/))' , 'L' , bc_velocityGrad ( : , : , j )
print '(a,/,3(3(f12.6,x)/))' , 'bc_stress' , bc_stress ( : , : , j )
print '(a,/,3(3(l,x)/))' , 'bc_mask' , bc_mask ( : , : , 1 , j )
print * , 'time' , bc_timeIncrement ( j )
print * , 'incs' , bc_steps ( j )
2010-06-10 20:21:10 +05:30
print * , ''
enddo
2010-06-25 17:01:05 +05:30
!read header of mesh file
a = 1_pInt
b = 1_pInt
c = 1_pInt
x = 1_pReal
y = 1_pReal
z = 1_pReal
2010-07-02 19:40:36 +05:30
gotResolution = . false .
gotDimension = . false .
gotHomogenization = . false .
2010-06-25 17:01:05 +05:30
path = getSolverJobName ( )
2010-07-01 20:50:06 +05:30
print * , 'JobName: ' , trim ( path )
2010-07-02 19:40:36 +05:30
if ( . not . IO_open_file ( unit , trim ( path ) / / InputFileExtension ) ) call IO_error ( 101 , ext_msg = path )
2010-06-25 17:01:05 +05:30
rewind ( unit )
do
2010-07-02 19:40:36 +05:30
read ( unit , '(a1024)' , END = 100 ) line
2010-06-25 17:01:05 +05:30
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
posMesh = IO_stringPos ( line , maxNchunksMesh )
select case ( IO_lc ( IO_StringValue ( line , posMesh , 1 ) ) )
case ( 'dimension' )
gotDimension = . true .
do i = 2 , 6 , 2
select case ( IO_lc ( IO_stringValue ( line , posMesh , i ) ) )
case ( 'x' )
x = IO_floatValue ( line , posMesh , i + 1 )
case ( 'y' )
y = IO_floatValue ( line , posMesh , i + 1 )
case ( 'z' )
z = IO_floatValue ( line , posMesh , i + 1 )
end select
enddo
case ( 'homogenization' )
gotHomogenization = . true .
homog = IO_intValue ( line , posMesh , 2 )
case ( 'resolution' )
gotResolution = . true .
do i = 2 , 6 , 2
select case ( IO_lc ( IO_stringValue ( line , posMesh , i ) ) )
case ( 'a' )
a = 2 ** IO_intValue ( line , posMesh , i + 1 )
case ( 'b' )
b = 2 ** IO_intValue ( line , posMesh , i + 1 )
case ( 'c' )
c = 2 ** IO_intValue ( line , posMesh , i + 1 )
end select
enddo
end select
if ( gotDimension . and . gotHomogenization . and . gotResolution ) exit
enddo
2010-07-01 20:50:06 +05:30
100 close ( unit )
2010-07-02 19:40:36 +05:30
print '(a,/,i3,i3,i3)' , 'resolution a b c' , a , b , c
print '(a,/,f6.2,f6.2,f6.2)' , 'dimension x y z' , x , y , z
print * , 'homogenization' , homog
print * , ''
2010-06-08 15:38:15 +05:30
2010-06-10 21:02:06 +05:30
temperature = 30 0.0_pReal
2010-07-01 20:50:06 +05:30
2010-07-02 19:40:36 +05:30
!-------------------------
!begin RL
!-------------------------
allocate ( datafft ( 2 * a * b * c ) )
allocate ( workfft ( 3 , 3 , a , b , c ) )
allocate ( workfftim ( 3 , 3 , a , b , c ) )
allocate ( sg ( 3 , 3 , a , b , c ) )
allocate ( disgrad ( 3 , 3 , a , b , c ) )
allocate ( defgradold ( 3 , 3 , a , b , c ) )
error = 0.00001
itmax = 100
delt ( 1 ) = 1.
delt ( 2 ) = 1.
delt ( 3 ) = 1.
nn ( 1 ) = a
nn ( 2 ) = b
nn ( 3 ) = c
nn2 ( 1 ) = a
nn2 ( 2 ) = b
2010-07-01 20:50:06 +05:30
2010-07-02 19:40:36 +05:30
prodnn = nn ( 1 ) * nn ( 2 ) * nn ( 3 )
wgt = 1. / prodnn
2010-07-01 20:50:06 +05:30
! C_0 and S_0 CALCULATION
!!
!! PHILIP: FE_exec_elem?
!!
2010-07-02 19:40:36 +05:30
stress = 0. !should be initialized somewhere
dsde = 0. !should be initialized somewhere
c0 = 0. !stiffness of reference material
c066 = 0. ! other way of notating c0
2010-07-05 17:03:48 +05:30
!#
!# do ielem = 1, prodnn !call each element with identity (math_i3) to initialize with high stress
!#
do ielem = 1 , int ( prodnn ) !call each element with identity (math_i3) to initialize with high stress
!#
2010-07-02 19:40:36 +05:30
call CPFEM_general ( 3 , math_i3 , math_i3 , temperature , 0.0_pReal , ielem , 1_pInt , stress , dsde )
2010-07-02 22:45:53 +05:30
enddo
2010-07-05 17:03:48 +05:30
!#
!# do ielem = 1, prodnn !call each element with identity (math_i3) to initialize with high stress
!#
do ielem = 1 , int ( prodnn ) !call each element with identity (math_i3) to initialize with high stress
!#
2010-07-02 22:45:53 +05:30
call CPFEM_general ( 2 , math_i3 , math_i3 , temperature , 0.0_pReal , ielem , 1_pInt , stress , dsde )
2010-07-02 19:40:36 +05:30
c066 = c066 + dsde
c0 = c0 + math_Mandel66to3333 ( dsde )
2010-07-01 20:50:06 +05:30
enddo
2010-07-02 19:40:36 +05:30
c066 = c066 * wgt
c0 = c0 * wgt
2010-07-01 20:50:06 +05:30
call math_invert ( 6 , c066 , s066 , idum , errmatinv )
if ( errmatinv ) then
2010-07-02 19:40:36 +05:30
write ( * , * ) 'ERROR IN C0 INVERSION'
stop
2010-07-01 20:50:06 +05:30
endif
2010-07-02 19:40:36 +05:30
s0 = math_Mandel66to3333 ( s066 )
2010-07-01 20:50:06 +05:30
! INITIALIZATION BEFORE STARTING WITH LOADINGS
2010-07-02 22:45:53 +05:30
disgrad = 0.0_pReal
disgradmacro = 0.0_pReal
2010-07-01 20:50:06 +05:30
2010-07-02 19:40:36 +05:30
do jload = 1 , N !Loop over loadcases defined in the loadcase file
udot ( : , : ) = bc_velocityGrad ( : , : , jload )
scauchy ( : , : ) = bc_stress ( : , : , jload )
iudot = 0
iscau = 0
do i = 1 , 3 !convert information about rb's from bc_mask in corresponding arrays
do j = 1 , 3
if ( bc_mask ( i , j , 1 , jload ) ) iudot ( i , j ) = 1
if ( bc_mask ( i , j , 2 , jload ) ) iscau ( i , j ) = 1
enddo
enddo
tdot = bc_timeIncrement ( jload ) / bc_steps ( jload )
! evm = dvm*tdot ?
do imicro = 1 , bc_steps ( jload ) ! loop oper steps defined in input file for current loadcase
write ( * , * ) '***************************************************'
write ( * , * ) 'STEP = ' , imicro
2010-07-01 20:50:06 +05:30
! INITIALIZATION BEFORE NEW TIME STEP
2010-07-02 19:40:36 +05:30
disgradmacro = disgradmacro + udot * tdot !update macroscopic displacementgradient
ddisgradmacro = 0.
ielem = 0
do k = 1 , c !loop over FPs
do j = 1 , b
do i = 1 , a
ielem = ielem + 1
2010-07-02 22:45:53 +05:30
defgradold ( : , : , i , j , k ) = math_I3 ( : , : ) + disgrad ( : , : , i , j , k ) ! wind forward
disgrad ( : , : , i , j , k ) = disgradmacro ( : , : ) ! no fluctuations as guess
call CPFEM_general ( 3 , defgradold ( : , : , i , j , k ) , math_i3 ( : , : ) + disgrad ( : , : , i , j , k ) , &
temperature , 0.0_pReal , ielem , 1_pInt , &
stress , dsde )
2010-07-01 20:50:06 +05:30
! sg(:,:,i,j,k)=math_Mandel6to33(stress) ?
2010-07-02 19:40:36 +05:30
enddo
enddo
enddo
2010-07-02 22:45:53 +05:30
ielem = 0
2010-07-02 19:40:36 +05:30
do k = 1 , c !loop over FPs
do j = 1 , b
do i = 1 , a
ielem = ielem + 1
2010-07-02 22:45:53 +05:30
call CPFEM_general ( 1 , defgradold ( : , : , i , j , k ) , math_i3 ( : , : ) + disgrad ( : , : , i , j , k ) , &
temperature , 0.0_pReal , ielem , 1_pInt , &
stress , dsde )
2010-07-02 19:40:36 +05:30
sg ( : , : , i , j , k ) = math_Mandel6to33 ( stress )
enddo
enddo
enddo
2010-07-01 20:50:06 +05:30
2010-07-02 22:45:53 +05:30
ddisgradmacroacum = 0.0_pReal
2010-07-02 19:40:36 +05:30
2010-07-02 22:45:53 +05:30
iter = 0_pInt
2010-07-05 17:03:48 +05:30
!# erre = 2.*error
2010-07-02 19:40:36 +05:30
errs = 2. * error
2010-07-05 17:03:48 +05:30
!# do while(iter <= itmax.and.(errs > error .or. erre > error))
do while ( iter < = itmax . and . errs > error )
2010-07-02 19:40:36 +05:30
iter = iter + 1
2010-07-02 22:45:53 +05:30
write ( * , * ) 'ITER = ' , iter
2010-07-02 19:40:36 +05:30
write ( * , * ) 'DIRECT FFT OF STRESS FIELD'
do ii = 1 , 3
do jj = 1 , 3
k1 = 0
do k = 1 , c
2010-07-01 20:50:06 +05:30
! write(*,'(1H+,a,i2,2(a,i4))')
! # 'STRESS - COMPONENT',ii,jj,' - Z = ',k,' OUT OF ',npts
2010-07-02 19:40:36 +05:30
do j = 1 , b
do i = 1 , a
k1 = k1 + 1
datafft ( k1 ) = sg ( ii , jj , i , j , k )
k1 = k1 + 1
datafft ( k1 ) = 0.
enddo
enddo
enddo
2010-07-02 22:45:53 +05:30
if ( c > 1 ) then
2010-07-02 19:40:36 +05:30
call fourn ( datafft , nn , 3 , 1 )
else
call fourn ( datafft , nn2 , 2 , 1 )
endif
k1 = 0
do k = 1 , c
do j = 1 , b
do i = 1 , a
k1 = k1 + 1
workfft ( ii , jj , i , j , k ) = datafft ( k1 )
k1 = k1 + 1
workfftim ( ii , jj , i , j , k ) = datafft ( k1 )
enddo
enddo
enddo
enddo
enddo
write ( * , * ) 'CALCULATING G^pq,ij : SG^ij ...'
do kzz = 1 , c
do kyy = 1 , b
do kxx = 1 , a
if ( kxx . le . a / 2 ) kx = kxx - 1
if ( kxx . gt . a / 2 ) kx = kxx - a - 1
if ( kyy . le . b / 2 ) ky = kyy - 1
if ( kyy . gt . b / 2 ) ky = kyy - b - 1
if ( kzz . le . c / 2 ) kz = kzz - 1
if ( kzz . gt . c / 2 ) kz = kzz - c - 1
xk ( 1 ) = kx / ( delt ( 1 ) * nn ( 1 ) )
xk ( 2 ) = ky / ( delt ( 2 ) * nn ( 2 ) )
if ( c . gt . 1 ) then
xk ( 3 ) = kz / ( delt ( 3 ) * nn ( 3 ) )
else
xk ( 3 ) = 0.
endif
2010-07-02 22:45:53 +05:30
! if (any(xk /= 0.0_pReal) then
! xknorm = 1.0_pReal/(xk(1)**2+xk(2)**2+xk(3)**2)
! forall (i=1:3,j=1:3) xkdyad(i,j) = xknorm * xk(i)*xk(j) ! the dyad is always used anfd could speed up things by using element-wise multiplication plus summation of array
! endif
2010-07-02 19:40:36 +05:30
xknorm = sqrt ( xk ( 1 ) ** 2 + xk ( 2 ) ** 2 + xk ( 3 ) ** 2 )
2010-07-02 22:45:53 +05:30
if ( xknorm / = 0.0_pReal ) then
2010-07-02 19:40:36 +05:30
do i = 1 , 3
2010-07-05 17:03:48 +05:30
!! xk2(i) = xk(i)/(xknorm*xknorm*2.*pi)
2010-07-02 19:40:36 +05:30
xk ( i ) = xk ( i ) / xknorm
enddo
endif
2010-07-02 22:45:53 +05:30
! forall(i=1:3,k=1:3) aux33(i,k) = sum(c0(i,:,k,:)*xkdyad(:,:)
! this is probably equiv with below quad looping
aux33 = 0.0_pReal
2010-07-02 19:40:36 +05:30
do i = 1 , 3
do k = 1 , 3
do j = 1 , 3
do l = 1 , 3
aux33 ( i , k ) = aux33 ( i , k ) + c0 ( i , j , k , l ) * xk ( j ) * xk ( l )
enddo
enddo
enddo
enddo
aux33 = math_inv3x3 ( aux33 )
2010-07-02 22:45:53 +05:30
2010-07-01 20:50:06 +05:30
! call minv(aux33,3,det,minv1,minv2)
! if(det.eq.0) then
! write(*,*) kx,ky,kz,' --> SINGULAR SYSTEM'
! stop
! pause
! endif
2010-07-02 22:45:53 +05:30
! forall (p=1:3,q=1:3,i=1:3,j=1:3) g1(p,q,i,j) = -aux33(p,i)*xk(q)*xk(j)
! could substitute below quad loop
2010-07-02 19:40:36 +05:30
do p = 1 , 3
do q = 1 , 3
do i = 1 , 3
do j = 1 , 3
g1 ( p , q , i , j ) = - aux33 ( p , i ) * xk ( q ) * xk ( j )
enddo
enddo
enddo
enddo
do i = 1 , 3
do j = 1 , 3
ddisgrad ( i , j ) = 0.
ddisgradim ( i , j ) = 0.
2010-07-01 20:50:06 +05:30
! if(kx.eq.0.and.ky.eq.0.and.kz.eq.0.) goto 4
2010-07-02 22:45:53 +05:30
if ( . not . ( kx == 0. . and . ky == 0. . and . kz == 0. ) ) then
2010-07-02 19:40:36 +05:30
do k = 1 , 3
do l = 1 , 3
ddisgrad ( i , j ) = ddisgrad ( i , j ) + g1 ( i , j , k , l ) * workfft ( k , l , kxx , kyy , kzz )
ddisgradim ( i , j ) = ddisgradim ( i , j ) + g1 ( i , j , k , l ) * workfftim ( k , l , kxx , kyy , kzz )
enddo
enddo
endif
enddo
enddo
workfft ( : , : , kxx , kyy , kzz ) = ddisgrad ( : , : )
workfftim ( : , : , kxx , kyy , kzz ) = ddisgradim ( : , : )
2010-07-02 22:45:53 +05:30
2010-07-02 19:40:36 +05:30
enddo
enddo
2010-07-01 20:50:06 +05:30
enddo
2010-07-02 19:40:36 +05:30
write ( * , * ) 'INVERSE FFT TO GET DISPLACEMENT GRADIENT FIELD'
do m1 = 1 , 3
do n1 = 1 , 3
k1 = 0
do k = 1 , c
do j = 1 , b
do i = 1 , a
2010-07-02 22:45:53 +05:30
2010-07-02 19:40:36 +05:30
k1 = k1 + 1
datafft ( k1 ) = workfft ( m1 , n1 , i , j , k )
2010-07-02 22:45:53 +05:30
2010-07-02 19:40:36 +05:30
k1 = k1 + 1
datafft ( k1 ) = workfftim ( m1 , n1 , i , j , k )
enddo
enddo
enddo
2010-07-02 22:45:53 +05:30
if ( c > 1 ) then
2010-07-02 19:40:36 +05:30
call fourn ( datafft , nn , 3 , - 1 )
else
call fourn ( datafft , nn2 , 2 , - 1 )
endif
datafft = datafft * wgt
k1 = 0
do k = 1 , c
do j = 1 , b
do i = 1 , a
k1 = k1 + 1
disgrad ( m1 , n1 , i , j , k ) = disgrad ( m1 , n1 , i , j , k ) + ddisgradmacro ( m1 , n1 ) + datafft ( k1 )
k1 = k1 + 1
enddo
enddo
enddo
enddo
enddo
2010-07-01 20:50:06 +05:30
2010-07-02 19:40:36 +05:30
write ( * , * ) 'UPDATE STRESS FIELD'
2010-07-01 20:50:06 +05:30
!!! call evpal
2010-07-02 19:40:36 +05:30
ielem = 0
do k = 1 , c
do j = 1 , b
do i = 1 , a
2010-07-01 20:50:06 +05:30
2010-07-02 19:40:36 +05:30
ielem = ielem + 1
2010-07-02 22:45:53 +05:30
call CPFEM_general ( 3 , defgradold ( : , : , i , j , k ) , math_i3 ( : , : ) + disgrad ( : , : , i , j , k ) , &
temperature , 0.0_pReal , ielem , 1_pInt , &
stress , dsde )
2010-07-01 20:50:06 +05:30
2010-07-02 19:40:36 +05:30
enddo
enddo
enddo
2010-07-01 20:50:06 +05:30
2010-07-02 19:40:36 +05:30
ielem = 0
2010-07-02 22:45:53 +05:30
scauav = 0.
2010-07-05 17:03:48 +05:30
!#
errs = 0.
!#
2010-07-02 19:40:36 +05:30
do k = 1 , c
do j = 1 , b
do i = 1 , a
2010-07-01 20:50:06 +05:30
2010-07-02 19:40:36 +05:30
ielem = ielem + 1
2010-07-02 22:45:53 +05:30
call CPFEM_general ( 2 , defgradold ( : , : , i , j , k ) , math_i3 ( : , : ) + disgrad ( : , : , i , j , k ) , &
temperature , 0.0_pReal , ielem , 1_pInt , &
stress , dsde )
2010-07-05 17:03:48 +05:30
!#
!# sg(:,:,i,j,k) = math_Mandel6to33(stress)
!# scauav(:,:) = scauav(:,:)+sg(:,:,i,j,k) ! average stress
!#
aux33 = math_Mandel6to33 ( stress )
erraux = 0.
do ii = 1 , 3
do jj = 1 , 3
erraux = erraux + ( sg ( ii , jj , i , j , k ) - aux33 ( ii , jj ) ) ** 2
enddo
enddo
errs = errs + sqrt ( erraux )
sg ( : , : , i , j , k ) = aux33
scauav = scauav + aux33 ! average stress
!#
2010-07-02 19:40:36 +05:30
enddo
enddo
enddo
2010-07-01 20:50:06 +05:30
2010-07-02 22:45:53 +05:30
scauav = scauav * wgt ! final weighting
2010-07-05 17:03:48 +05:30
!#
errs = errs * wgt
scaunorm = 0.
do ii = 1 , 3
do jj = 1 , 3
scaunorm = scaunorm + scauav ( ii , jj ) ** 2
enddo
enddo
scaunorm = sqrt ( scaunorm )
errs = errs / scaunorm
!#
2010-07-01 20:50:06 +05:30
! MIXED BC
2010-07-02 19:40:36 +05:30
do i = 1 , 3
do j = 1 , 3
ddisgradmacro ( i , j ) = 0.
2010-07-02 22:45:53 +05:30
if ( iudot ( i , j ) == 0 ) then
! ddisgradmacro(i,j) = ddisgradmacro(i,j)+sum(s0(i,j,:,:)*iscau(:,:)*(scauchy(:,:)-scauav(:,:)))=
! could replace the k,l loop
2010-07-02 19:40:36 +05:30
do k = 1 , 3
do l = 1 , 3
ddisgradmacro ( i , j ) = ddisgradmacro ( i , j ) + s0 ( i , j , k , l ) * iscau ( k , l ) * ( scauchy ( k , l ) - scauav ( k , l ) )
enddo
enddo
endif
enddo
enddo
ddisgradmacroacum = ddisgradmacroacum + ddisgradmacro
! write(*,*) 'DDEFGRADMACRO(1,1),(2,2) = ',ddefgradmacro(1,1),ddefgradmacro(2,2)
!! svm = ?
! erre = erre/evm
! errs = errs/svm
write ( * , * ) 'STRESS FIELD ERROR = ' , errs
write ( * , * ) 'STRAIN FIELD ERROR = ' , erre
2010-07-01 20:50:06 +05:30
! write(21,101) iter,erre,errs,svm
!101 format(i3,4(1x,e10.4),10(1x,F7.4))
2010-07-02 19:40:36 +05:30
enddo ! WHILE ENDDO
2010-07-01 20:50:06 +05:30
2010-07-02 19:40:36 +05:30
disgradmacroactual = disgradmacro + ddisgradmacroacum
2010-07-01 20:50:06 +05:30
!! write(*,*) 'defgradmacro(1,1),defgradmacro(2,2),defgradmacro(3,3)'
!! write(*,*) defgradmacroactual(1,1),defgradmacroactual(2,2),defgradmacroactual(3,3)
!! write(*,*) 'defgradmacro(1,1)/defgradmacro(3,3)'
!! write(*,*) defgradmacroactual(1,1)/defgradmacroactual(3,3)
!! write(*,*) 'scauav(1,1),scauav(2,2),scauav(3,3)'
!! write(*,*) scauav(1,1),scauav(2,2),scauav(3,3)
2010-07-02 19:40:36 +05:30
enddo ! IMICRO ENDDO
enddo ! JLOAD ENDDO Ende looping over loadcases
2010-07-01 20:50:06 +05:30
!-------------------------
2010-07-02 19:40:36 +05:30
!end RL
2010-07-01 20:50:06 +05:30
!-------------------------
2010-06-10 20:21:10 +05:30
end program mpie_spectral
!********************************************************************
! quit subroutine to satisfy IO_error
!
!********************************************************************
2010-06-08 15:38:15 +05:30
subroutine quit ( id )
use prec
implicit none
integer ( pInt ) id
stop
end subroutine
2010-07-01 20:50:06 +05:30
2010-07-02 19:40:36 +05:30
!********************************************************************
! fourn subroutine (fourier transform)
2010-07-01 20:50:06 +05:30
! FROM NUMERICAL RECIPES IN F77 (FIXED FORMAT),
! CONVERTED INTO FREE FORMAT (RL @ MPIE, JUNE 2010)
2010-07-02 19:40:36 +05:30
!********************************************************************
subroutine fourn ( data , nn , ndim , isign )
2010-07-01 20:50:06 +05:30
INTEGER isign , ndim , nn ( ndim )
REAL data ( * )
! INTEGER i1,i2,i2rev,i3,i3rev,ibit,idim,ifp1,ifp2,ip1,ip2,ip3,k1,
! *k2,n,nprev,nrem,ntot
INTEGER i1 , i2 , i2rev , i3 , i3rev , ibit , idim , ifp1 , ifp2 , ip1 , ip2 , ip3 , k1 , k2 , n , nprev , nrem , ntot
REAL tempi , tempr
DOUBLE PRECISION theta , wi , wpi , wpr , wr , wtemp
2010-07-02 19:40:36 +05:30
ntot = 1
! do 11 idim = 1,ndim
do idim = 1 , ndim ! 11
ntot = ntot * nn ( idim )
2010-07-01 20:50:06 +05:30
!11 continue
enddo ! 11
!
2010-07-02 19:40:36 +05:30
nprev = 1
! do 18 idim = 1,ndim
do idim = 1 , ndim ! 18
n = nn ( idim )
nrem = ntot / ( n * nprev )
ip1 = 2 * nprev
ip2 = ip1 * n
ip3 = ip2 * nrem
i2rev = 1
! do 14 i2 = 1,ip2,ip1
do i2 = 1 , ip2 , ip1 ! 14
2010-07-01 20:50:06 +05:30
if ( i2 . lt . i2rev ) then
2010-07-02 19:40:36 +05:30
! do 13 i1 = i2,i2+ip1-2,2
! do 12 i3 = i1,ip3,ip2
do i1 = i2 , i2 + ip1 - 2 , 2 ! 13
do i3 = i1 , ip3 , ip2 ! 12
i3rev = i2rev + i3 - i2
tempr = data ( i3 )
tempi = data ( i3 + 1 )
data ( i3 ) = data ( i3rev )
data ( i3 + 1 ) = data ( i3rev + 1 )
data ( i3rev ) = tempr
data ( i3rev + 1 ) = tempi
2010-07-01 20:50:06 +05:30
!13 continue
enddo ! 12
enddo ! 13
endif
2010-07-02 19:40:36 +05:30
ibit = ip2 / 2
2010-07-01 20:50:06 +05:30
!1 if ((ibit.ge.ip1).and.(i2rev.gt.ibit)) then
do while ( ( ibit . ge . ip1 ) . and . ( i2rev . gt . ibit ) ) ! if 1
2010-07-02 19:40:36 +05:30
i2rev = i2rev - ibit
ibit = ibit / 2
2010-07-01 20:50:06 +05:30
! goto 1
! endif
enddo ! do while (if 1)
2010-07-02 19:40:36 +05:30
i2rev = i2rev + ibit
2010-07-01 20:50:06 +05:30
!14 continue
enddo ! 14
2010-07-02 19:40:36 +05:30
ifp1 = ip1
2010-07-01 20:50:06 +05:30
!2 if(ifp1.lt.ip2)then
do while ( ifp1 . lt . ip2 ) ! if 2
2010-07-02 19:40:36 +05:30
ifp2 = 2 * ifp1
theta = isign * 6.28318530717959d0 / ( ifp2 / ip1 )
wpr = - 2.d0 * sin ( 0.5d0 * theta ) ** 2
wpi = sin ( theta )
wr = 1.d0
wi = 0.d0
! do 17 i3 = 1,ifp1,ip1
! do 16 i1 = i3,i3+ip1-2,2
! do 15 i2 = i1,ip3,ifp2
do i3 = 1 , ifp1 , ip1 ! 17
do i1 = i3 , i3 + ip1 - 2 , 2 ! 16
do i2 = i1 , ip3 , ifp2 ! 15
k1 = i2
k2 = k1 + ifp1
tempr = sngl ( wr ) * data ( k2 ) - sngl ( wi ) * data ( k2 + 1 )
tempi = sngl ( wr ) * data ( k2 + 1 ) + sngl ( wi ) * data ( k2 )
data ( k2 ) = data ( k1 ) - tempr
data ( k2 + 1 ) = data ( k1 + 1 ) - tempi
data ( k1 ) = data ( k1 ) + tempr
data ( k1 + 1 ) = data ( k1 + 1 ) + tempi
2010-07-01 20:50:06 +05:30
!15 continue
!16 continue
enddo ! 15
enddo ! 16
2010-07-02 19:40:36 +05:30
wtemp = wr
wr = wr * wpr - wi * wpi + wr
wi = wi * wpr + wtemp * wpi + wi
2010-07-01 20:50:06 +05:30
!17 continue
enddo ! 17
2010-07-02 19:40:36 +05:30
ifp1 = ifp2
2010-07-01 20:50:06 +05:30
! goto 2
! endif
enddo ! do while (if 2)
2010-07-02 19:40:36 +05:30
nprev = n * nprev
2010-07-01 20:50:06 +05:30
!18 continue
enddo ! 18
return
END subroutine