F77 --> F90 polishing/condensation as far as possible.
next step is change of Re/Im FFT to Re only...
This commit is contained in:
parent
4d110126da
commit
49926d5d66
|
@ -241,6 +241,8 @@ program mpie_spectral
|
||||||
use IO
|
use IO
|
||||||
use math
|
use math
|
||||||
use CPFEM, only: CPFEM_general
|
use CPFEM, only: CPFEM_general
|
||||||
|
use FEsolving, only: FEsolving_execElem, FEsolving_execIP
|
||||||
|
use debug
|
||||||
|
|
||||||
implicit none
|
implicit none
|
||||||
|
|
||||||
|
@ -285,16 +287,16 @@ program mpie_spectral
|
||||||
real(pReal), dimension(3,3) :: disgradmacro, disgradmacroactual
|
real(pReal), dimension(3,3) :: disgradmacro, disgradmacroactual
|
||||||
real(pReal), dimension(3,3) :: ddisgradmacro, ddisgradmacroacum, ddisgrad, ddisgradim
|
real(pReal), dimension(3,3) :: ddisgradmacro, ddisgradmacroacum, ddisgrad, ddisgradim
|
||||||
real(pReal), dimension(3,3) :: defgrad0, defgrad
|
real(pReal), dimension(3,3) :: defgrad0, defgrad
|
||||||
real(pReal), dimension(3,3) :: udot, scauchy, scauav, aux33, xkdyad
|
real(pReal), dimension(3,3) :: udot, scauchy, scauav, aux33, xkdyad, xknormdyad
|
||||||
|
|
||||||
!integer(pInt), dimension(2) :: nn2 m.diehl
|
!integer(pInt), dimension(2) :: nn2 m.diehl
|
||||||
|
|
||||||
real(pReal), dimension(3) :: delt,xk
|
real(pReal), dimension(3) :: delt,xk
|
||||||
real(pReal), dimension(6) :: aux6
|
real(pReal), dimension(6) :: aux6
|
||||||
|
|
||||||
integer(pInt) prodnn,itmax, jload, ielem, ii, jj, k1, kxx, kyy, kzz, kx, ky, kz, idum, iter, imicro, m1, n1, p, q
|
integer(pInt) prodnn,itmax, jload, ielem, ii, jj, k1, kxx, kyy, kzz, kx, ky, kz, idum, iter, imicro, m1, n1, p, q, ione
|
||||||
|
|
||||||
real(pReal) wgt,error,tdot,erre,errs,evm,svm,det,xknorm, erraux, scaunorm
|
real(pReal) wgt,error,timestep,erre,errs,evm,svm,det,xknorm, erraux, scaunorm
|
||||||
logical errmatinv
|
logical errmatinv
|
||||||
|
|
||||||
|
|
||||||
|
@ -380,14 +382,15 @@ program mpie_spectral
|
||||||
|
|
||||||
! consistency checks
|
! consistency checks
|
||||||
do j = 1,N
|
do j = 1,N
|
||||||
if (any(bc_mask(:,:,1,j) == bc_mask(:,:,2,j))) &
|
! if (any(bc_mask(:,:,1,j) == bc_mask(:,:,2,j))) & ! don't enforce consistency to allow values as initial guess
|
||||||
call IO_error(47,j) ! bc_mask consistency
|
! call IO_error(47,j) ! bc_mask consistency
|
||||||
if (any(math_transpose3x3(bc_stress(:,:,j)) + bc_stress(:,:,j) /= 2.0_pReal * bc_stress(:,:,j))) &
|
if (any(math_transpose3x3(bc_stress(:,:,j)) + bc_stress(:,:,j) /= 2.0_pReal * bc_stress(:,:,j))) &
|
||||||
call IO_error(48,j) ! bc_stress symmetry
|
call IO_error(48,j) ! bc_stress symmetry
|
||||||
|
|
||||||
print '(a,/,3(3(f12.6,x)/))','L',bc_velocityGrad(:,:,j)
|
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(f12.6,x)/))','bc_stress',bc_stress(:,:,j)
|
||||||
print '(a,/,3(3(l,x)/))','bc_mask',bc_mask(:,:,1,j)
|
print '(a,/,3(3(l,x)/))','bc_mask for velocitygrad',bc_mask(:,:,1,j)
|
||||||
|
print '(a,/,3(3(l,x)/))','bc_mask for stress',bc_mask(:,:,2,j)
|
||||||
print *,'time',bc_timeIncrement(j)
|
print *,'time',bc_timeIncrement(j)
|
||||||
print *,'incs',bc_steps(j)
|
print *,'incs',bc_steps(j)
|
||||||
print *, ''
|
print *, ''
|
||||||
|
@ -453,14 +456,25 @@ program mpie_spectral
|
||||||
|
|
||||||
|
|
||||||
allocate (datafft(2*resolution(1)*resolution(2)*resolution(3)))
|
allocate (datafft(2*resolution(1)*resolution(2)*resolution(3)))
|
||||||
allocate (workfft(3,3,resolution(1),resolution(2),resolution(3)))
|
!allocate (datafft(resolution(1)*resolution(2)*resolution(3))) !for real fft m.diehl
|
||||||
allocate (workfftim(3,3,resolution(1),resolution(2),resolution(3)))
|
|
||||||
allocate (sg(3,3,resolution(1),resolution(2),resolution(3)))
|
allocate (workfft(resolution(1),resolution(2),resolution(3),3,3))
|
||||||
|
allocate (workfftim(resolution(1),resolution(2),resolution(3),3,3)) ! probably not needed for real fft m.diehl
|
||||||
|
allocate (sg(resolution(1),resolution(2),resolution(3),3,3))
|
||||||
allocate (disgrad(3,3,resolution(1),resolution(2),resolution(3)))
|
allocate (disgrad(3,3,resolution(1),resolution(2),resolution(3)))
|
||||||
allocate (defgradold(3,3,resolution(1),resolution(2),resolution(3)))
|
allocate (defgradold(3,3,resolution(1),resolution(2),resolution(3)))
|
||||||
|
|
||||||
|
open(56,file='str_str.out',status='unknown')
|
||||||
|
|
||||||
|
write(56,'(t2,a,t14,a,t26,a,t38,a,t50,a,t62,a,t74,a,t86,a,t98,a,&
|
||||||
|
t110,a,t122,a,t134,a,t146,a,t158,a,t170,a)')&
|
||||||
|
'U1,1','U2,2','U3,3', 'U2,3','U3,1','U1,2','U3,2','U1,3','U2,1',&
|
||||||
|
'S11','S22','S33','S23','S31','S12'
|
||||||
|
|
||||||
|
ione=1
|
||||||
|
|
||||||
error = 0.00001
|
error = 0.00001
|
||||||
itmax = 100
|
itmax = 50
|
||||||
|
|
||||||
delt(1) = 1.
|
delt(1) = 1.
|
||||||
delt(2) = 1.
|
delt(2) = 1.
|
||||||
|
@ -504,27 +518,35 @@ program mpie_spectral
|
||||||
|
|
||||||
do jload = 1,N !Loop over loadcases defined in the loadcase file
|
do jload = 1,N !Loop over loadcases defined in the loadcase file
|
||||||
udot(:,:) = bc_velocityGrad(:,:,jload)
|
udot(:,:) = bc_velocityGrad(:,:,jload)
|
||||||
|
! udot(1,1)=-.35 ! temporary, to solve problem with bc mask
|
||||||
|
! udot(2,2)=-.35 not needed any more
|
||||||
scauchy(:,:) = bc_stress(:,:,jload)
|
scauchy(:,:) = bc_stress(:,:,jload)
|
||||||
iudot = 0
|
iudot = 0
|
||||||
iscau = 0
|
iscau = 0
|
||||||
|
|
||||||
do i = 1, 3 !convert information about rb's from bc_mask in corresponding arrays
|
do i = 1, 3 !convert information about rb's from bc_mask in corresponding arrays
|
||||||
do j = 1, 3
|
do j = 1, 3
|
||||||
if(bc_mask(i,j,1,jload)) iudot(i,j) = 1
|
! if(bc_mask(i,j,1,jload)) iudot(i,j) = 1
|
||||||
if(bc_mask(i,j,2,jload)) iscau(i,j) = 1
|
! if(bc_mask(i,j,2,jload)) iscau(i,j) = 1!
|
||||||
enddo
|
if(bc_mask(i,j,2,jload)) then
|
||||||
enddo
|
iscau(i,j) = 1
|
||||||
|
iudot(i,j) = 0
|
||||||
|
else
|
||||||
|
iscau(i,j) = 0
|
||||||
|
iudot(i,j) = 1
|
||||||
|
endif
|
||||||
|
enddo; enddo
|
||||||
|
|
||||||
|
timestep = bc_timeIncrement(jload)/bc_steps(jload)
|
||||||
|
|
||||||
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
|
do imicro = 1, bc_steps(jload) ! loop oper steps defined in input file for current loadcase
|
||||||
write(*,*) '***************************************************'
|
write(*,*) '***************************************************'
|
||||||
write(*,*) 'STEP = ',imicro
|
write(*,*) 'STEP = ',imicro
|
||||||
|
|
||||||
! INITIALIZATION BEFORE NEW TIME STEP
|
! INITIALIZATION BEFORE NEW TIME STEP
|
||||||
disgradmacro = disgradmacro+udot*tdot !update macroscopic displacementgradient
|
disgradmacro = disgradmacro+udot*timestep !update macroscopic displacementgradient
|
||||||
ddisgradmacro = 0.
|
ddisgradmacro = 0._pReal
|
||||||
ielem=0
|
ielem = 0_pInt
|
||||||
|
|
||||||
do k = 1, resolution(3) !loop over FPs
|
do k = 1, resolution(3) !loop over FPs
|
||||||
do j = 1, resolution(2)
|
do j = 1, resolution(2)
|
||||||
|
@ -533,25 +555,25 @@ program mpie_spectral
|
||||||
defgradold(:,:,i,j,k) = math_I3(:,:) + disgrad(:,:,i,j,k) ! wind forward
|
defgradold(:,:,i,j,k) = math_I3(:,:) + disgrad(:,:,i,j,k) ! wind forward
|
||||||
disgrad(:,:,i,j,k) = disgradmacro(:,:) ! no fluctuations as guess
|
disgrad(:,:,i,j,k) = disgradmacro(:,:) ! no fluctuations as guess
|
||||||
call CPFEM_general(3,defgradold(:,:,i,j,k),math_I3(:,:)+disgrad(:,:,i,j,k),&
|
call CPFEM_general(3,defgradold(:,:,i,j,k),math_I3(:,:)+disgrad(:,:,i,j,k),&
|
||||||
temperature,0.0_pReal,ielem,1_pInt,&
|
temperature,timestep,ielem,1_pInt,&
|
||||||
stress,dsde, pstress, dPdF)
|
stress,dsde, pstress, dPdF)
|
||||||
enddo
|
enddo; enddo; enddo
|
||||||
enddo
|
|
||||||
enddo
|
|
||||||
|
|
||||||
ielem = 0
|
ielem = 0
|
||||||
|
call debug_reset()
|
||||||
|
|
||||||
do k = 1, resolution(3) !loop over FPs
|
do k = 1, resolution(3) !loop over FPs
|
||||||
do j = 1, resolution(2)
|
do j = 1, resolution(2)
|
||||||
do i = 1, resolution(1)
|
do i = 1, resolution(1)
|
||||||
ielem = ielem+1
|
ielem = ielem+1
|
||||||
call CPFEM_general(1,defgradold(:,:,i,j,k),math_i3(:,:)+disgrad(:,:,i,j,k),&
|
call CPFEM_general(min(2,ielem),& ! first element gets calcMode 1, others 2 (saves winding forward effort)
|
||||||
temperature,0.0_pReal,ielem,1_pInt,&
|
defgradold(:,:,i,j,k),math_i3(:,:)+disgrad(:,:,i,j,k),&
|
||||||
|
temperature,timestep,ielem,1_pInt,&
|
||||||
stress,dsde, pstress, dPdF)
|
stress,dsde, pstress, dPdF)
|
||||||
sg(:,:,i,j,k) = math_Mandel6to33(stress)
|
sg(i,j,k,:,:) = math_Mandel6to33(stress)
|
||||||
enddo
|
enddo; enddo; enddo
|
||||||
enddo
|
|
||||||
enddo
|
call debug_info()
|
||||||
|
|
||||||
ddisgradmacroacum = 0.0_pReal
|
ddisgradmacroacum = 0.0_pReal
|
||||||
|
|
||||||
|
@ -560,283 +582,200 @@ program mpie_spectral
|
||||||
errs = 2.*error
|
errs = 2.*error
|
||||||
|
|
||||||
!# do while(iter <= itmax.and.(errs > error .or. erre > error))
|
!# do while(iter <= itmax.and.(errs > error .or. erre > error))
|
||||||
do while(iter <= itmax .and. errs > error)
|
do while(iter < itmax .and. errs > error)
|
||||||
iter = iter+1
|
iter = iter+1
|
||||||
write(*,*) 'ITER = ',iter
|
write(*,*) 'ITER = ',iter
|
||||||
write(*,*) 'DIRECT FFT OF STRESS FIELD'
|
write(*,*) 'DIRECT FFT OF STRESS FIELD'
|
||||||
do ii = 1,3
|
do ii = 1,3
|
||||||
do jj = 1,3
|
do jj = 1,3
|
||||||
k1 = 0
|
datafft = 0._pReal
|
||||||
do k = 1, resolution(3)
|
datafft(1:2*prodnn:2) = reshape(sg(:,:,:,ii,jj),(/prodnn/))
|
||||||
! write(*,'(1H+,a,i2,2(a,i4))')
|
|
||||||
! # 'STRESS - COMPONENT',ii,jj,' - Z = ',k,' OUT OF ',npts
|
|
||||||
do j = 1, resolution(2)
|
|
||||||
do i = 1, resolution(1)
|
|
||||||
k1 = k1+1
|
|
||||||
datafft(k1) = sg(ii,jj,i,j,k)
|
|
||||||
k1 = k1+1
|
|
||||||
datafft(k1) = 0.
|
|
||||||
enddo
|
|
||||||
enddo
|
|
||||||
enddo
|
|
||||||
if(resolution(3) > 1) then
|
if(resolution(3) > 1) then
|
||||||
call fourn(datafft,resolution,3,1)
|
call fourn(datafft,resolution,3,1)
|
||||||
else
|
else
|
||||||
call fourn(datafft,resolution(1:2),2,1) !nn2 replaced by resolution(1:2)
|
call fourn(datafft,resolution(1:2),2,1)
|
||||||
endif
|
endif
|
||||||
k1 = 0
|
|
||||||
do k = 1, resolution(3)
|
workfft(:,:,:,ii,jj) = reshape(datafft(1:2*prodnn:2),resolution)
|
||||||
do j = 1, resolution(2)
|
workfftim(:,:,:,ii,jj) = reshape(datafft(2:2*prodnn:2),resolution)
|
||||||
do i = 1, resolution(1)
|
|
||||||
k1 = k1+1
|
enddo; enddo
|
||||||
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 ...'
|
write(*,*) 'CALCULATING G^pq,ij : SG^ij ...'
|
||||||
|
|
||||||
do kzz = 1, resolution(3)
|
do kzz = 1, resolution(3)
|
||||||
|
kz = kzz-1
|
||||||
|
if(kzz > resolution(3)/2) kz = kz-resolution(3)
|
||||||
|
if(resolution(3) > 1) then
|
||||||
|
xk(3) = kz/(delt(3)*resolution(3))
|
||||||
|
else
|
||||||
|
xk(3) = 0.
|
||||||
|
endif
|
||||||
do kyy = 1, resolution(2)
|
do kyy = 1, resolution(2)
|
||||||
|
ky = kyy-1
|
||||||
|
if(kyy > resolution(2)/2) ky = ky-resolution(2)
|
||||||
|
xk(2) = ky/(delt(2)*resolution(2))
|
||||||
do kxx = 1, resolution(1)
|
do kxx = 1, resolution(1)
|
||||||
if(kxx <= resolution(1)/2) kx = kxx-1
|
kx = kxx-1
|
||||||
if(kxx > resolution(1)/2) kx = kxx-resolution(1)-1
|
if(kxx > resolution(1)/2) kx = kx-resolution(1)
|
||||||
|
|
||||||
if(kyy <= resolution(2)/2) ky = kyy-1
|
|
||||||
if(kyy > resolution(2)/2) ky = kyy-resolution(2)-1
|
|
||||||
|
|
||||||
if(kzz <= resolution(3)/2) kz = kzz-1
|
|
||||||
if(kzz > resolution(3)/2) kz = kzz-resolution(3)-1
|
|
||||||
|
|
||||||
xk(1) = kx/(delt(1)*resolution(1))
|
xk(1) = kx/(delt(1)*resolution(1))
|
||||||
xk(2) = ky/(delt(2)*resolution(2))
|
|
||||||
|
|
||||||
if(resolution(3) > 1) then
|
|
||||||
xk(3) = kz/(delt(3)*resolution(3))
|
|
||||||
else
|
|
||||||
xk(3) = 0.
|
|
||||||
endif
|
|
||||||
|
|
||||||
|
forall (i=1:3,j=1:3) xkdyad(i,j) = xk(i)*xk(j) ! the dyad is always used and could speed up things by using element-wise multiplication plus summation of array
|
||||||
if (any(xk /= 0.0_pReal)) then
|
if (any(xk /= 0.0_pReal)) then
|
||||||
xknorm = 1.0_pReal/(xk(1)**2+xk(2)**2+xk(3)**2)
|
xknormdyad = xkdyad * 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 and could speed up things by using element-wise multiplication plus summation of array
|
else
|
||||||
|
xknormdyad = xkdyad
|
||||||
endif
|
endif
|
||||||
forall(i=1:3,k=1:3) aux33(i,k) = sum(c0(i,:,k,:)*xkdyad(:,:))
|
forall(i=1:3,k=1:3) aux33(i,k) = sum(c0(i,:,k,:)*xknormdyad(:,:))
|
||||||
|
|
||||||
|
|
||||||
!!!!!!!!! replace that by the stuff above, except from small rounding errors get the same results m.diehl
|
|
||||||
! xknorm = sqrt(xk(1)**2+xk(2)**2+xk(3)**2)
|
|
||||||
!
|
|
||||||
! if (xknorm /= 0.0_pReal) then
|
|
||||||
! do i = 1,3
|
|
||||||
! ! xk2(i) = xk(i)/(xknorm*xknorm*2.*pi)
|
|
||||||
! xk(i) = xk(i)/xknorm
|
|
||||||
! enddo
|
|
||||||
! endif
|
|
||||||
!
|
|
||||||
! aux33 = 0.0_pReal
|
|
||||||
! 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)
|
aux33 = math_inv3x3(aux33)
|
||||||
|
forall (p=1:3,q=1:3,i=1:3,j=1:3) g1(p,q,i,j) = -aux33(p,i)*xkdyad(q,j)
|
||||||
|
|
||||||
! call minv(aux33,3,det,minv1,minv2)
|
ddisgrad = 0._pReal
|
||||||
! if(det.eq.0) then
|
ddisgradim = 0._pReal
|
||||||
! write(*,*) kx,ky,kz,' --> SINGULAR SYSTEM'
|
|
||||||
! stop
|
|
||||||
! pause
|
|
||||||
! endif
|
|
||||||
|
|
||||||
! 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
|
|
||||||
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 i = 1,3
|
||||||
do j = 1,3
|
do j = 1,3
|
||||||
ddisgrad(i,j) = 0.
|
|
||||||
ddisgradim(i,j) = 0.
|
|
||||||
! if(kx.eq.0.and.ky.eq.0.and.kz.eq.0.) goto 4
|
|
||||||
if(kx /= 0 .or. ky /= 0 .or. kz /= 0) then
|
if(kx /= 0 .or. ky /= 0 .or. kz /= 0) then
|
||||||
do k = 1,3
|
ddisgrad(i,j) = ddisgrad(i,j) + sum(g1(i,j,:,:)*workfft(kxx,kyy,kzz,:,:))
|
||||||
do l = 1,3
|
ddisgradim(i,j) = ddisgradim(i,j) + sum(g1(i,j,:,:)*workfftim(kxx,kyy,kzz,:,:))
|
||||||
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
|
endif
|
||||||
enddo
|
enddo; enddo
|
||||||
enddo
|
|
||||||
|
|
||||||
workfft(:,:,kxx,kyy,kzz) = ddisgrad(:,:)
|
workfft(kxx,kyy,kzz,:,:) = ddisgrad(:,:)
|
||||||
workfftim(:,:,kxx,kyy,kzz) = ddisgradim(:,:)
|
workfftim(kxx,kyy,kzz,:,:) = ddisgradim(:,:)
|
||||||
|
|
||||||
enddo
|
enddo; enddo; enddo
|
||||||
enddo
|
|
||||||
enddo
|
|
||||||
|
|
||||||
write(*,*) 'INVERSE FFT TO GET DISPLACEMENT GRADIENT FIELD'
|
write(*,*) 'INVERSE FFT TO GET DISPLACEMENT GRADIENT FIELD'
|
||||||
|
|
||||||
do m1 = 1,3
|
do ii = 1,3
|
||||||
do n1 = 1,3
|
do jj = 1,3
|
||||||
|
|
||||||
k1 = 0
|
datafft(1:2*prodnn:2) = reshape(workfft(:,:,:,ii,jj),(/prodnn/))
|
||||||
|
datafft(2:2*prodnn:2) = reshape(workfftim(:,:,:,ii,jj),(/prodnn/))
|
||||||
do k = 1, resolution(3)
|
|
||||||
do j = 1, resolution(2)
|
|
||||||
do i = 1, resolution(1)
|
|
||||||
|
|
||||||
k1 = k1+1
|
|
||||||
datafft(k1) = workfft(m1,n1,i,j,k)
|
|
||||||
|
|
||||||
k1 = k1+1
|
|
||||||
datafft(k1) = workfftim(m1,n1,i,j,k)
|
|
||||||
enddo
|
|
||||||
enddo
|
|
||||||
enddo
|
|
||||||
|
|
||||||
if(resolution(3) > 1) then ! distinguish 2D and 3D case
|
if(resolution(3) > 1) then ! distinguish 2D and 3D case
|
||||||
call fourn(datafft,resolution,3,-1)
|
call fourn(datafft,resolution,3,-1)
|
||||||
else
|
else
|
||||||
call fourn(datafft,resolution(1:2),2,-1) !nn2 replaced by resolution(1:2)
|
call fourn(datafft,resolution(1:2),2,-1)
|
||||||
endif
|
endif
|
||||||
datafft = datafft*wgt
|
|
||||||
|
|
||||||
k1 = 0
|
disgrad(ii,jj,:,:,:) = disgrad(ii,jj,:,:,:) + ddisgradmacro(ii,jj)
|
||||||
do k = 1, resolution(3)
|
disgrad(ii,jj,:,:,:) = disgrad(ii,jj,:,:,:) + reshape(datafft(1:2*prodnn:2)*wgt,resolution)
|
||||||
do j = 1, resolution(2)
|
|
||||||
do i = 1, resolution(1)
|
|
||||||
|
|
||||||
k1 = k1+1
|
enddo; enddo
|
||||||
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
|
|
||||||
|
|
||||||
write(*,*) 'UPDATE STRESS FIELD'
|
write(*,*) 'UPDATE STRESS FIELD'
|
||||||
!!! call evpal
|
|
||||||
ielem = 0
|
ielem = 0_pInt
|
||||||
do k = 1, resolution(3)
|
do k = 1, resolution(3)
|
||||||
do j = 1, resolution(2)
|
do j = 1, resolution(2)
|
||||||
do i = 1, resolution(1)
|
do i = 1, resolution(1)
|
||||||
|
|
||||||
ielem = ielem+1
|
ielem = ielem+1
|
||||||
call CPFEM_general(3,defgradold(:,:,i,j,k),math_i3(:,:)+disgrad(:,:,i,j,k),&
|
call CPFEM_general(3,defgradold(:,:,i,j,k),math_i3(:,:)+disgrad(:,:,i,j,k),&
|
||||||
temperature,0.0_pReal,ielem,1_pInt,&
|
temperature,timestep,ielem,1_pInt,&
|
||||||
stress,dsde, pstress, dPdF)
|
stress,dsde, pstress, dPdF)
|
||||||
|
|
||||||
enddo
|
enddo; enddo; enddo
|
||||||
enddo
|
|
||||||
enddo
|
ielem = 0_pInt
|
||||||
|
scauav = 0._pReal
|
||||||
|
errs = 0._pReal
|
||||||
|
|
||||||
|
call debug_reset()
|
||||||
|
|
||||||
ielem = 0
|
|
||||||
scauav = 0.
|
|
||||||
!#
|
|
||||||
errs=0.
|
|
||||||
!#
|
|
||||||
do k = 1, resolution(3)
|
do k = 1, resolution(3)
|
||||||
do j = 1, resolution(2)
|
do j = 1, resolution(2)
|
||||||
do i = 1, resolution(1)
|
do i = 1, resolution(1)
|
||||||
|
|
||||||
ielem = ielem+1
|
ielem = ielem+1
|
||||||
call CPFEM_general(2,defgradold(:,:,i,j,k),math_i3(:,:)+disgrad(:,:,i,j,k),&
|
call CPFEM_general(2,defgradold(:,:,i,j,k),math_i3(:,:)+disgrad(:,:,i,j,k),&
|
||||||
temperature,0.0_pReal,ielem,1_pInt,&
|
temperature,timestep,ielem,1_pInt,&
|
||||||
stress,dsde, pstress, dPdF)
|
stress,dsde, pstress, dPdF)
|
||||||
!#
|
|
||||||
!# sg(:,:,i,j,k) = math_Mandel6to33(stress)
|
|
||||||
!# scauav(:,:) = scauav(:,:)+sg(:,:,i,j,k) ! average stress
|
|
||||||
!#
|
|
||||||
aux33 = math_Mandel6to33(stress)
|
|
||||||
|
|
||||||
erraux=0.
|
aux33 = math_Mandel6to33(stress)
|
||||||
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
|
errs = errs + sqrt(sum((sg(i,j,k,:,:)-aux33(:,:))**2))
|
||||||
scauav = scauav+aux33 ! average stress
|
|
||||||
!#
|
sg(i,j,k,:,:) = aux33
|
||||||
enddo
|
scauav = scauav + aux33 ! average stress
|
||||||
enddo
|
|
||||||
enddo
|
enddo; enddo; enddo
|
||||||
|
|
||||||
|
call debug_info()
|
||||||
|
|
||||||
|
errs = errs/sqrt(sum(scauav(:,:)**2))
|
||||||
|
|
||||||
scauav = scauav*wgt ! final weighting
|
scauav = scauav*wgt ! final weighting
|
||||||
!#
|
|
||||||
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
|
|
||||||
!#
|
|
||||||
|
|
||||||
! MIXED BC
|
! MIXED BC
|
||||||
|
|
||||||
|
|
||||||
|
ddisgradmacro = 0._pReal
|
||||||
do i = 1,3
|
do i = 1,3
|
||||||
do j = 1,3
|
do j = 1,3
|
||||||
ddisgradmacro(i,j) = 0.
|
if(iudot(i,j) == 0) &
|
||||||
if(iudot(i,j) == 0) then
|
ddisgradmacro(i,j) = sum(s0(i,j,:,:)*iscau(:,:)*(scauchy(:,:)-scauav(:,:)))
|
||||||
! ddisgradmacro(i,j) = ddisgradmacro(i,j)+sum(s0(i,j,:,:)*iscau(:,:)*(scauchy(:,:)-scauav(:,:)))=
|
enddo; enddo
|
||||||
! could replace the k,l loop
|
|
||||||
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
|
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(*,*) 'STRESS FIELD ERROR = ',errs
|
||||||
write(*,*) 'STRAIN FIELD ERROR = ',erre
|
write(*,*) 'STRAIN FIELD ERROR = ',erre
|
||||||
! write(21,101) iter,erre,errs,svm
|
! write(21,101) iter,erre,errs,svm
|
||||||
!101 format(i3,4(1x,e10.4),10(1x,F7.4))
|
!101 format(i3,4(1x,e10.4),10(1x,F7.4))
|
||||||
|
|
||||||
enddo ! WHILE ENDDO
|
enddo ! convergence iteration
|
||||||
|
|
||||||
disgradmacroactual = disgradmacro+ddisgradmacroacum
|
disgradmacroactual = disgradmacro+ddisgradmacroacum
|
||||||
|
|
||||||
!! write(*,*) 'defgradmacro(1,1),defgradmacro(2,2),defgradmacro(3,3)'
|
write(*,*) 'U1,1,U2,2,U3,3'
|
||||||
!! write(*,*) defgradmacroactual(1,1),defgradmacroactual(2,2),defgradmacroactual(3,3)
|
write(*,*) disgradmacroactual(1,1),disgradmacroactual(2,2),disgradmacroactual(3,3)
|
||||||
!! write(*,*) 'defgradmacro(1,1)/defgradmacro(3,3)'
|
write(*,*) 'U1,1/U3,3'
|
||||||
!! write(*,*) defgradmacroactual(1,1)/defgradmacroactual(3,3)
|
write(*,*) disgradmacroactual(1,1)/disgradmacroactual(3,3)
|
||||||
!! write(*,*) 'scauav(1,1),scauav(2,2),scauav(3,3)'
|
write(*,*) 'S11,S22,S33'
|
||||||
!! write(*,*) scauav(1,1),scauav(2,2),scauav(3,3)
|
write(*,*) scauav(1,1),scauav(2,2),scauav(3,3)
|
||||||
|
|
||||||
enddo ! IMICRO ENDDO
|
write(56,'(15(e11.4,1x))') disgradmacroactual(1,1),disgradmacroactual(2,2),disgradmacroactual(3,3), &
|
||||||
enddo ! JLOAD ENDDO Ende looping over loadcases
|
disgradmacroactual(2,3),disgradmacroactual(3,1),disgradmacroactual(1,2), &
|
||||||
|
disgradmacroactual(3,2),disgradmacroactual(1,3),disgradmacroactual(2,1), &
|
||||||
|
scauav(1,1),scauav(2,2),scauav(3,3), &
|
||||||
|
scauav(2,3),scauav(3,1),scauav(1,2)
|
||||||
|
|
||||||
|
IF(IMICRO.EQ.1.OR.IMICRO.EQ.40) THEN
|
||||||
|
|
||||||
|
IF(IMICRO.EQ.1) THEN
|
||||||
|
open(91,file='fields1.out',status='unknown')
|
||||||
|
ELSE IF(IMICRO.EQ.40) THEN
|
||||||
|
open(91,file='fields40.out',status='unknown')
|
||||||
|
ENDIF
|
||||||
|
|
||||||
|
write(91,*) delt
|
||||||
|
|
||||||
|
write(91,*) ' x y z ngr ph Ui,j ... Sij ...'
|
||||||
|
|
||||||
|
do k=1,resolution(3)
|
||||||
|
do j=1,resolution(2)
|
||||||
|
do i=1,resolution(1)
|
||||||
|
|
||||||
|
!! write(91,'(5i5,18(e11.3,1x))') i,j,k,jgrain(i,j,k),jphase(i,j,k),
|
||||||
|
write(91,'(5i5,18(e11.3,1x))') i,j,k,ione,ione, &
|
||||||
|
disgrad(1,1,i,j,k),disgrad(2,2,i,j,k),disgrad(3,3,i,j,k), &
|
||||||
|
disgrad(2,3,i,j,k),disgrad(3,1,i,j,k),disgrad(1,2,i,j,k), &
|
||||||
|
disgrad(3,2,i,j,k),disgrad(1,3,i,j,k),disgrad(2,1,i,j,k), &
|
||||||
|
sg(i,j,k,1,1),sg(i,j,k,2,2),sg(i,j,k,3,3), &
|
||||||
|
sg(i,j,k,2,3),sg(i,j,k,3,1),sg(i,j,k,1,2)
|
||||||
|
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
|
|
||||||
|
close(91)
|
||||||
|
|
||||||
|
ENDIF
|
||||||
|
|
||||||
|
enddo ! time stepping through increment
|
||||||
|
enddo ! loadcases
|
||||||
|
|
||||||
end program mpie_spectral
|
end program mpie_spectral
|
||||||
|
|
||||||
|
@ -863,6 +802,7 @@ end subroutine
|
||||||
subroutine fourn(data,nn,ndim,isign)
|
subroutine fourn(data,nn,ndim,isign)
|
||||||
|
|
||||||
use prec, only: pInt,pReal
|
use prec, only: pInt,pReal
|
||||||
|
use math, only: pi
|
||||||
implicit none
|
implicit none
|
||||||
|
|
||||||
integer(pInt) isign,ndim,nn(ndim)
|
integer(pInt) isign,ndim,nn(ndim)
|
||||||
|
@ -909,18 +849,18 @@ subroutine fourn(data,nn,ndim,isign)
|
||||||
|
|
||||||
do while (ifp1.lt.ip2)
|
do while (ifp1.lt.ip2)
|
||||||
ifp2 = 2*ifp1
|
ifp2 = 2*ifp1
|
||||||
theta = isign*6.28318530717959d0/(ifp2/ip1)
|
theta = isign*2_pReal*pi/(ifp2/ip1)
|
||||||
wpr = -2.d0*sin(0.5d0*theta)**2
|
wpr = -2_pReal*sin(0.5_pReal*theta)**2
|
||||||
wpi = sin(theta)
|
wpi = sin(theta)
|
||||||
wr = 1.d0
|
wr = 1_pReal
|
||||||
wi = 0.d0
|
wi = 0_pReal
|
||||||
do i3 = 1,ifp1,ip1 ! 17
|
do i3 = 1,ifp1,ip1 ! 17
|
||||||
do i1 = i3,i3+ip1-2,2 ! 16
|
do i1 = i3,i3+ip1-2,2 ! 16
|
||||||
do i2 = i1,ip3,ifp2 ! 15
|
do i2 = i1,ip3,ifp2 ! 15
|
||||||
k1 = i2
|
k1 = i2
|
||||||
k2 = k1+ifp1
|
k2 = k1+ifp1
|
||||||
tempr = sngl(wr)*data(k2)-sngl(wi)*data(k2+1)
|
tempr = wr*data(k2)-wi*data(k2+1)
|
||||||
tempi = sngl(wr)*data(k2+1)+sngl(wi)*data(k2)
|
tempi = wr*data(k2+1)+wi*data(k2)
|
||||||
data(k2) = data(k1)-tempr
|
data(k2) = data(k1)-tempr
|
||||||
|
|
||||||
data(k2+1) = data(k1+1)-tempi
|
data(k2+1) = data(k1+1)-tempi
|
||||||
|
|
Loading…
Reference in New Issue