changed declaration of two variables from real to int, spectral method is now working.

started to commend Ricardos code, layoutet loops and removed some redundant variables.
until now, no error calculation is done. at the moment calculations are in an infinite loop
This commit is contained in:
Martin Diehl 2010-07-02 14:10:36 +00:00
parent 661bb97800
commit d6ba9d54b6
1 changed files with 395 additions and 489 deletions

View File

@ -269,18 +269,11 @@ program mpie_spectral
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
!-------------------------
!begin RL
!-------------------------
! begin change by m.diehl
!integer(pInt), parameter :: npts1 = 32
!integer(pInt), parameter :: npts2 = 32
!integer(pInt), parameter :: npts3 = 32
integer(pInt) npts1, npts2, npts3 ! will be allocated later
! end change by m.diehl
real(pReal), dimension(:), allocatable :: datafft
real(pReal), dimension(:,:,:,:,:), allocatable :: workfft,workfftim,sg,disgrad,defgradold
@ -291,13 +284,15 @@ program mpie_spectral
real(pReal), dimension(3,3) :: defgrad0,defgrad
real(pReal), dimension(3,3) :: udot,scauchy,scauav,aux33
real(pReal), dimension(3) :: delt,nn,xk,xk2
real(pReal), dimension(2) :: nn2
integer(pInt), dimension(3) :: nn
integer(pInt), dimension(2) :: nn2
real(pReal), dimension(3) :: delt,xk,xk2
real(pReal), dimension(6) :: aux6
real(pReal), dimension(3,3,3,3) :: c0,s0,g1
real(pReal), dimension(6,6) :: c066,s066
integer(pInt) itmax,nsteps,jload,ielem,ii,jj,k1,kxx,kyy,kzz,kx,ky,kz,idum,iter,imicro,m1,n1,p,q
integer(pInt) itmax, jload, ielem, ii, jj, k1, kxx, kyy, kzz, kx, ky, kz, idum, iter, imicro, m1, n1, p, q
real(pReal) prodnn,wgt,error,tdot,erre,errs,evm,svm,det,xknorm
@ -461,21 +456,16 @@ program mpie_spectral
print *, ''
temperature = 300.0_pReal
!-------------------------
!RL
!-------------------------
!begin change m.diehl
npts1=a
npts2=b
npts3=c
!end change m.diehl
allocate (datafft(2*npts1*npts2*npts3))
allocate (workfft(3,3,npts1,npts2,npts3))
allocate (workfftim(3,3,npts1,npts2,npts3))
allocate (sg(3,3,npts1,npts2,npts3))
allocate (disgrad(3,3,npts1,npts2,npts3))
allocate (defgradold(3,3,npts1,npts2,npts3))
!-------------------------
!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
@ -484,46 +474,36 @@ program mpie_spectral
delt(2) = 1.
delt(3) = 1.
nn(1)=npts1
nn(2)=npts2
nn(3)=npts3
nn(1) = a
nn(2) = b
nn(3) = c
nn2(1)=npts1
nn2(2)=npts2
nn2(1) = a
nn2(2) = b
prodnn = nn(1)*nn(2)*nn(3)
wgt = 1./prodnn
! C_0 and S_0 CALCULATION
!!
!! PHILIP: FE_exec_elem?
!!
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
ielem=0
c0=0.
c066=0.
do k=1,npts3
do j=1,npts2
do i=1,npts1
ielem=ielem+1
do ielem = 1, prodnn !call each element with identity (math_i3) to initialize with high stress
call CPFEM_general(3,math_i3,math_i3,temperature,0.0_pReal,ielem,1_pInt,stress,dsde)
c066 = c066+dsde
c0 = c0+math_Mandel66to3333(dsde)
enddo
enddo
enddo
c066 = c066*wgt
c0 = c0*wgt
call math_invert(6,c066,s066,idum,errmatinv)
call math_invert(6,c066,s066,idum,errmatinv)
if(errmatinv) then
write(*,*) 'ERROR IN C0 INVERSION'
stop
@ -535,84 +515,64 @@ program mpie_spectral
disgradmacro = 0.
do k=1,npts3
do j=1,npts2
do i=1,npts1
do k = 1, c
do j = 1, b
do i = 1, a
defgradold(:,:,i,j,k) = math_i3(:,:)
enddo
enddo
enddo
do jload = 1,N
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
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
nsteps=bc_steps(jload)
tdot = bc_timeIncrement(jload)/bc_steps(jload)
! evm = dvm*tdot ?
do imicro=1,nsteps
do imicro = 1, bc_steps(jload) ! loop oper steps defined in input file for current loadcase
write(*,*) '***************************************************'
write(*,*) 'STEP = ',imicro
! write(21,*) 'STEP = ',imicro
! INITIALIZATION BEFORE NEW TIME STEP
disgradmacro=disgradmacro+udot*tdot
disgradmacro = disgradmacro+udot*tdot !update macroscopic displacementgradient
ddisgradmacro = 0.
ielem=0
do k=1,npts3
do j=1,npts2
do i=1,npts1
do k = 1, c !loop over FPs
do j = 1, b
do i = 1, a
ielem = ielem+1
disgrad(:,:,i,j,k)=disgradmacro(:,:)
defgrad0(:,:)=defgradold(:,:,i,j,k)
defgrad(:,:)=math_i3(:,:)+disgrad(:,:,i,j,k)
defgrad(:,:)=math_i3(:,:)+disgrad(:,:,i,j,k) ! calculate new defgrad
call CPFEM_general(3,defgrad0,defgrad,temperature,0.0_pReal,ielem,1_pInt,stress,dsde)
! sg(:,:,i,j,k)=math_Mandel6to33(stress) ?
enddo
enddo
enddo
ielem = 0
do k=1,npts3
do j=1,npts2
do i=1,npts1
do k = 1, c !loop over FPs
do j = 1, b
do i = 1, a
ielem = ielem+1
disgrad(:,:,i,j,k) = disgradmacro(:,:)
defgrad0(:,:) = defgradold(:,:,i,j,k)
defgrad(:,:) = math_i3(:,:)+disgrad(:,:,i,j,k)
call CPFEM_general(1,defgrad0,defgrad,temperature,0.0_pReal,ielem,1_pInt,stress,dsde)
sg(:,:,i,j,k) = math_Mandel6to33(stress)
enddo
enddo
enddo
@ -623,79 +583,60 @@ program mpie_spectral
erre = 2.*error
errs = 2.*error
do while(iter.lt.itmax.and.(errs.gt.error.or.erre.gt.error))
do while(iter <= itmax.and.(errs.gt.error.or.erre.gt.error))
iter = iter+1
write(*,*)'ITER = ',iter
write(*,*) 'DIRECT FFT OF STRESS FIELD'
do ii = 1,3
do jj = 1,3
k1 = 0
do k=1,npts3
do k = 1,c
! write(*,'(1H+,a,i2,2(a,i4))')
! # 'STRESS - COMPONENT',ii,jj,' - Z = ',k,' OUT OF ',npts
do j=1,npts2
do i=1,npts1
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
if(npts3.gt.1) then
if(c.gt.1) then
call fourn(datafft,nn,3,1)
else
call fourn(datafft,nn2,2,1)
endif
k1 = 0
do k=1,npts3
do j=1,npts2
do i=1,npts1
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
do kzz=1,npts3
do kyy=1,npts2
do kxx=1,npts1
if(kyy.le.b/2) ky = kyy-1
if(kyy.gt.b/2) ky = kyy-b-1
if(kxx.le.npts1/2) kx=kxx-1
if(kxx.gt.npts1/2) kx=kxx-npts1-1
if(kyy.le.npts2/2) ky=kyy-1
if(kyy.gt.npts2/2) ky=kyy-npts2-1
if(kzz.le.npts3/2) kz=kzz-1
if(kzz.gt.npts3/2) kz=kzz-npts3-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(npts3.gt.1) then
if(c.gt.1) then
xk(3) = kz/(delt(3)*nn(3))
else
xk(3) = 0.
@ -720,7 +661,6 @@ program mpie_spectral
enddo
enddo
enddo
aux33 = math_inv3x3(aux33)
! call minv(aux33,3,det,minv1,minv2)
@ -742,24 +682,17 @@ program mpie_spectral
do i = 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(.not.(kx.eq.0.and.ky.eq.0.and.kz.eq.0.)) then
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
@ -777,37 +710,34 @@ program mpie_spectral
k1 = 0
do k=1,npts3
do j=1,npts2
do i=1,npts1
do k = 1,c
do j = 1,b
do i = 1,a
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(npts3.gt.1) then
if(c.gt.1) then
call fourn(datafft,nn,3,-1)
else
call fourn(datafft,nn2,2,-1)
endif
datafft = datafft*wgt
k1 = 0
do k=1,npts3
do j=1,npts2
do i=1,npts1
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
@ -816,14 +746,11 @@ program mpie_spectral
enddo
write(*,*) 'UPDATE STRESS FIELD'
!!! call evpal
ielem = 0
do k=1,npts3
do j=1,npts2
do i=1,npts1
do k = 1,c
do j = 1,b
do i = 1,a
ielem = ielem+1
@ -833,16 +760,15 @@ program mpie_spectral
call CPFEM_general(3,defgrad0,defgrad,temperature,0.0_pReal,ielem,1_pInt,stress,dsde)
sg(:,:,i,j,k) = math_Mandel6to33(stress)
enddo
enddo
enddo
ielem = 0
do k=1,npts3
do j=1,npts2
do i=1,npts1
do k = 1,c
do j = 1,b
do i = 1,a
ielem = ielem+1
@ -852,61 +778,43 @@ program mpie_spectral
call CPFEM_general(2,defgrad0,defgrad,temperature,0.0_pReal,ielem,1_pInt,stress,dsde)
sg(:,:,i,j,k) = math_Mandel6to33(stress)
enddo
enddo
enddo
!!! call get_smacro
! AVERAGE STRESS
scauav = 0.
do k=1,npts3
do j=1,npts2
do i=1,npts1
do k = 1,c
do j = 1,b
do i = 1,a
scauav(:,:) = scauav(:,:)+sg(:,:,i,j,k)*wgt
enddo
enddo
enddo
! MIXED BC
do i = 1,3
do j = 1,3
ddisgradmacro(i,j) = 0.
if(iudot(i,j).eq.0) then
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
! write(21,101) iter,erre,errs,svm
!101 format(i3,4(1x,e10.4),10(1x,F7.4))
@ -921,24 +829,19 @@ program mpie_spectral
!! write(*,*) 'scauav(1,1),scauav(2,2),scauav(3,3)'
!! write(*,*) scauav(1,1),scauav(2,2),scauav(3,3)
do k=1,npts3
do j=1,npts2
do i=1,npts1
do k = 1,c
do j = 1,b
do i = 1,a
defgradold(:,:,i,j,k) = math_i3(:,:)+disgrad(:,:,i,j,k)
enddo
enddo
enddo
enddo ! IMICRO ENDDO
enddo ! JLOAD ENDDO
enddo ! JLOAD ENDDO Ende looping over loadcases
!-------------------------
!RL
!end RL
!-------------------------
end program mpie_spectral
@ -955,11 +858,15 @@ subroutine quit(id)
stop
end subroutine
SUBROUTINE fourn(data,nn,ndim,isign)
!
!********************************************************************
! fourn subroutine (fourier transform)
! FROM NUMERICAL RECIPES IN F77 (FIXED FORMAT),
! CONVERTED INTO FREE FORMAT (RL @ MPIE, JUNE 2010)
!
!********************************************************************
subroutine fourn(data,nn,ndim,isign)
INTEGER isign,ndim,nn(ndim)
REAL data(*)
! INTEGER i1,i2,i2rev,i3,i3rev,ibit,idim,ifp1,ifp2,ip1,ip2,ip3,k1,
@ -970,7 +877,6 @@ end subroutine
ntot = 1
! do 11 idim = 1,ndim
do idim = 1,ndim ! 11
!
ntot = ntot*nn(idim)
!11 continue
enddo ! 11
@ -999,7 +905,7 @@ end subroutine
data(i3+1) = data(i3rev+1)
data(i3rev) = tempr
data(i3rev+1) = tempi
!12 continue
!13 continue
enddo ! 12
enddo ! 13