Added random component

This commit is contained in:
Luc Hantcherli 2007-03-28 17:48:59 +00:00
parent 3ba46ac7d7
commit 11bb8faa3f
1 changed files with 142 additions and 117 deletions

View File

@ -1,3 +1,6 @@
! QUESTION fileunit 1 may run into trouble?
!************************************
!* Module: CONSTITUTIVE *
!************************************
@ -159,8 +162,7 @@ data constitutive_sd(:,12,3)/ 1, 1, 0/ ; data constitutive_sn(:,12,3)/ 1,-1, 1/
!* Slip-slip interactions matrices
!* (defined for the moment as crystal structure property and not as material property)
!* (may be changed in the future)
real(pReal), dimension(constitutive_MaxMaxNslipOfStructure,constitutive_MaxMaxNslipOfStructure,constitutive_MaxCrystalStructure) ::&
constitutive_HardeningMatrix
real(pReal), dimension(constitutive_MaxMaxNslipOfStructure,constitutive_MaxMaxNslipOfStructure,constitutive_MaxCrystalStructure) :: constitutive_HardeningMatrix
real(pReal), parameter :: constitutive_LatentHardening=1.4_pReal
!*************************************
@ -201,6 +203,7 @@ character(len=80), dimension(:), allocatable :: texture_symmetry
integer(pInt), dimension(:) , allocatable :: texture_Ngrains
integer(pInt), dimension(:) , allocatable :: texture_NGauss
integer(pInt),dimension(:) , allocatable :: texture_NFiber
integer(pInt),dimension(:) , allocatable :: texture_NRandom
real(pReal), dimension(:,:,:) , allocatable :: texture_Gauss
real(pReal), dimension(:,:,:) , allocatable :: texture_Fiber
real(pReal), dimension(:,:,:) , allocatable :: constitutive_phi1
@ -282,8 +285,7 @@ do l=1,3
constitutive_Sslip(i,j,k,l)=constitutive_sd(i,k,l)*constitutive_sn(j,k,l)
endforall
!* Normalization of Schmid matrix
invNorm=dsqrt(1.0_pReal/((constitutive_sn(1,k,l)**2+constitutive_sn(2,k,l)**2+constitutive_sn(3,k,l)**2)* &
(constitutive_sd(1,k,l)**2+constitutive_sd(2,k,l)**2+constitutive_sd(3,k,l)**2)))
invNorm=dsqrt(1.0_pReal/((constitutive_sn(1,k,l)**2+constitutive_sn(2,k,l)**2+constitutive_sn(3,k,l)**2)*(constitutive_sd(1,k,l)**2+constitutive_sd(2,k,l)**2+constitutive_sd(3,k,l)**2)))
constitutive_Sslip(:,:,k,l)=constitutive_Sslip(:,:,k,l)*invNorm
!* Vectorization of normalized Schmid matrix
!* according MARC component order 11,22,33,12,23,13
@ -367,7 +369,7 @@ integer(pInt), dimension(3) :: positions
count=0
part=''
do while(.true.)
do
read(file,'(a80)',END=100) line
positions=IO_stringPos(line,1)
tag=IO_lc(IO_stringValue(line,positions,1))
@ -383,7 +385,7 @@ enddo
end subroutine
subroutine constitutive_CountGaussAndFiber(file,count,part)
character(len=80) function constitutive_assignNGaussAndFiber(file)
!*********************************************************************
!*********************************************************************
use prec, only: pInt
@ -391,30 +393,33 @@ use IO, only: IO_stringPos,IO_stringValue,IO_lc
implicit none
!* Definition of variables
character(len=80) line,tag,part
integer(pInt) file,count,pos
character(len=80) line,tag
integer(pInt) file,section,pos
integer(pInt), dimension(3) :: positions
part=''
constitutive_assignNGaussAndFiber=''
section = 0_pInt
do while(.true.)
do
read(file,'(a80)',END=100) line
positions=IO_stringPos(line,1)
tag=IO_lc(IO_stringValue(line,positions,1))
if (tag(1:1)=='<'.AND.tag(len_trim(tag):len_trim(tag))=='>') then
part=tag(2:len_trim(tag)-1)
constitutive_assignNGaussAndFiber=tag(2:len_trim(tag)-1)
exit
elseif (tag(1:1)=='[') then
count=count+1
elseif (tag(2:len_trim(tag)-1)=='gauss') then
texture_NGauss(count)=texture_NGauss(count)+1
elseif (tag(2:len_trim(tag)-1)=='fiber') then
texture_NFiber(count)=texture_NFiber(count)+1
section=section+1
texture_NGauss(section) = 0_pInt
texture_NFiber(section) = 0_pInt
elseif (tag=='(gauss)') then
texture_NGauss(section)=texture_NGauss(section)+1
elseif (tag=='(fiber)') then
texture_NFiber(section)=texture_NFiber(section)+1
endif
enddo
100 return
end subroutine
end function
character(len=80) function constitutive_Parse_UnknownPart(file)
@ -436,7 +441,7 @@ integer(pInt), dimension(1+2*maxNchunks) :: positions
constitutive_parse_unknownPart=''
do while(.true.)
do
read(file,'(a80)',END=100) line
positions=IO_stringPos(line,maxNchunks)
tag=IO_lc(IO_stringValue(line,positions,1))
@ -572,9 +577,9 @@ do while(.true.)
case('phi2')
texture_Gauss(3,gaussCount,section)=IO_floatValue(line,positions,i+1)*inRad
case('scatter')
texture_Gauss(5,gaussCount,section)=IO_floatValue(line,positions,i+1)
texture_Gauss(4,gaussCount,section)=IO_floatValue(line,positions,i+1)
case('fraction')
texture_Gauss(6,gaussCount,section)=IO_floatValue(line,positions,i+1)
texture_Gauss(5,gaussCount,section)=IO_floatValue(line,positions,i+1)
end select
enddo
case ('(fiber)')
@ -616,16 +621,17 @@ subroutine constitutive_Parse_MatTexDat(filename)
!* - filename : name of input file *
!*********************************************************************
use prec, only: pReal,pInt
use IO
use IO, only: IO_error
!use math, only: math_Mandel3333to66, math_Voigt66to3333
implicit none
!* Definition of variables
character(len=*) filename
character(len=80) part,formerPart
integer(pInt) sectionCount,dummy,i,j,k
integer(pInt) sectionCount,i,j,k
!* First reading: number of materials and textures
!* determine material_maxN and texture_maxN
!* determine material_maxN and texture_maxN from last respective parts
open(1,FILE=filename,ACTION='READ',STATUS='OLD',ERR=100)
part = '_dummy_'
do while (part/='')
@ -638,28 +644,7 @@ do while (part/='')
texture_maxN = sectionCount
end select
enddo
close(1)
!* Arrays allocation
allocate(texture_NGauss(texture_maxN)) ; texture_NGauss=0_pInt
allocate(texture_NFiber(texture_maxN)) ; texture_NFiber=0_pInt
!* Second reading: number of Gauss and Fiber
!* determine material_maxN and texture_maxN
open(1,FILE=filename,ACTION='READ',STATUS='OLD',ERR=100)
part = '_dummy_'
sectionCount = 0
do while (part/='')
select case (part)
case ('textures')
call constitutive_CountGaussAndFiber(1,sectionCount,part)
case default
call constitutive_CountSections(1,dummy,part)
end select
enddo
close(1)
!* Arrays allocation
texture_maxNGauss=maxval(texture_NGauss)
texture_maxNFiber=maxval(texture_NFiber)
allocate(material_CrystalStructure(material_maxN)) ; material_CrystalStructure=0_pInt
allocate(material_Nslip(material_maxN)) ; material_Nslip=0_pInt
allocate(material_C11(material_maxN)) ; material_C11=0.0_pReal
@ -677,11 +662,29 @@ allocate(material_w0(material_maxN)) ; material_w0=0.0_pRea
allocate(texture_ODFfile(texture_maxN)) ; texture_ODFfile=''
allocate(texture_Ngrains(texture_maxN)) ; texture_Ngrains=0_pInt
allocate(texture_symmetry(texture_maxN)) ; texture_symmetry=''
allocate(texture_Gauss(6,texture_maxNGauss,texture_maxN)) ; texture_Gauss=0.0_pReal
allocate(texture_Fiber(6,texture_maxNGauss,texture_maxN)) ; texture_Fiber=0.0_pReal
allocate(texture_NGauss(texture_maxN)) ; texture_NGauss=0_pInt
allocate(texture_NFiber(texture_maxN)) ; texture_NFiber=0_pInt
allocate(texture_NRandom(texture_maxN)) ; texture_NRandom=0_pInt
!* Second reading: number of Gauss and Fiber
rewind(1)
part = '_dummy_'
do while (part/='')
select case (part)
case ('textures')
part = constitutive_assignNGaussAndFiber(1)
case default
part = constitutive_Parse_UnknownPart(1)
end select
enddo
!* Arrays allocation
texture_maxNGauss=maxval(texture_NGauss)
texture_maxNFiber=maxval(texture_NFiber)
allocate(texture_Gauss(5,texture_maxNGauss,texture_maxN)) ; texture_Gauss=0.0_pReal
allocate(texture_Fiber(6,texture_maxNFiber,texture_maxN)) ; texture_Fiber=0.0_pReal
!* Third reading: materials and textures are stored
open(1,FILE=filename,ACTION='READ',STATUS='OLD',ERR=100)
rewind(1)
part='_dummy_'
do while (part/='')
select case (part)
@ -699,7 +702,7 @@ close(1)
!* Construction of the elasticity matrices
do i=1,material_maxN
select case (material_CrystalStructure(i))
case(1:2)
case(1:2) ! cubic(s)
do k=1,3
do j=1,3
material_Cslip_66(k,j,i)=material_C12(i)
@ -707,7 +710,7 @@ do i=1,material_maxN
material_Cslip_66(k,k,i)=material_C11(i)
material_Cslip_66(k+3,k+3,i)=material_C44(i)
enddo
case(3)
case(3) ! hcp
material_Cslip_66(1,1,i)=material_C11(i)
material_Cslip_66(2,2,i)=material_C11(i)
material_Cslip_66(3,3,i)=material_C33(i)
@ -721,15 +724,14 @@ do i=1,material_maxN
material_Cslip_66(5,5,i)=material_C44(i)
material_Cslip_66(6,6,i)=0.5_pReal*(material_C11(i)-material_C12(i))
end select
!* If MARC
!constitutive_Cslip_66(:,:,i)=TransformCForMarc(constitutive_Cslip_66(:,:,i))
! material_Cslip_66(:,:,i) = math_Mandel3333to66(math_Voigt66to3333(material_Cslip_66(:,:,i)))
enddo
! MISSING some consistency checks may be..?
! if ODFfile present then set NGauss NFiber =0
return
100 call IO_error(110) ! corrupt matarials_textures file
100 call IO_error(110) ! corrupt materials_textures file
end subroutine
@ -743,26 +745,40 @@ use mesh, only: mesh_NcpElems,FE_Nips,mesh_maxNips,mesh_element
implicit none
!* Definition of variables
integer(pInt) i,j,k,l
integer(pInt) mul
integer(pInt) i,j,k,l,volfrac,Ngrains
integer(pInt) matID,texID,multiplicity
!* Check for random components
do i=1,texture_maxN
if (texture_ODFfile(i)=='') then
volfrac=sum(texture_gauss(5,:,i))+sum(texture_fiber(6,:,i))
if (volfrac<1.0_pReal) then
texture_NRandom(i)=1
endif
endif
enddo
!* Multiplicity of orientations per texture
!* Construction of optimized constitutive_Ngrains
allocate(constitutive_Ngrains(mesh_maxNips,mesh_NcpElems))
constitutive_Ngrains=0_pInt
allocate(constitutive_Ngrains(mesh_maxNips,mesh_NcpElems)) ; constitutive_Ngrains=0_pInt
do i=1,mesh_NcpElems
do j=1,FE_Nips(mesh_element(2,i))
mul=texture_Ngrains(mesh_element(4,i))/(texture_NGauss(mesh_element(4,i))+texture_NFiber(mesh_element(4,i)))
if (texture_Ngrains(mesh_element(4,i))==(mul*(texture_NGauss(mesh_element(4,i))+texture_NFiber(mesh_element(4,i))))) then
constitutive_Ngrains(j,i)=texture_Ngrains(mesh_element(4,i))
else
constitutive_Ngrains(j,i)=mul*(texture_NGauss(mesh_element(4,i))+texture_NFiber(mesh_element(4,i)))
endif
constitutive_maxNgrains=maxval(constitutive_Ngrains)
enddo
matID=mesh_element(3,i)
texID=mesh_element(4,i)
if (texture_ODFfile(texID)/='') then
constitutive_Ngrains = texture_Ngrains(texID)
else
multiplicity=texture_Ngrains(texID)/(texture_NGauss(texID)+texture_NFiber(texID)+texture_NRandom(texID))
do j=1,FE_Nips(mesh_element(2,i))
if (texture_Ngrains(texID)==multiplicity*(texture_NGauss(texID)+texture_NFiber(texID)+texture_NRandom(texID))) then
constitutive_Ngrains(j,i)=texture_Ngrains(texID)
else
constitutive_Ngrains(j,i)=multiplicity*(texture_NGauss(texID)+texture_NFiber(texID)+texture_NRandom(texID))
endif
enddo
endif
enddo
!* Allocate arrays
constitutive_maxNgrains=maxval(constitutive_Ngrains)
allocate(constitutive_matID(constitutive_maxNgrains,mesh_maxNips,mesh_NcpElems))
constitutive_matID=0_pInt
allocate(constitutive_texID(constitutive_maxNgrains,mesh_maxNips,mesh_NcpElems))
@ -788,45 +804,63 @@ constitutive_state_new=0.0_pReal
!* Results
allocate(constitutive_Nresults(constitutive_maxNgrains,mesh_maxNips,mesh_NcpElems))
constitutive_Nresults=0_pInt
!* Assignment
do i=1,mesh_NcpElems
do j=1,FE_Nips(mesh_element(2,i))
mul=texture_Ngrains(mesh_element(4,i))/(texture_NGauss(mesh_element(4,i))+texture_NFiber(mesh_element(4,i)))
!* Gauss component
do k=1,mul*texture_NGauss(mesh_element(4,i)),mul
do l=k,k+mul
constitutive_matID(l,j,i)=mesh_element(3,i)
constitutive_texID(l,j,i)=mesh_element(4,i)
constitutive_MatVolFrac(l,j,i)=1.0_pReal
constitutive_TexVolFrac(l,j,i)=texture_Gauss(6,k,mesh_element(4,i))/mul
constitutive_phi1(l,j,i)=texture_Gauss(1,k,mesh_element(4,i))
constitutive_phi(l,j,i)=texture_Gauss(2,k,mesh_element(4,i))
constitutive_phi2(l,j,i)=texture_Gauss(3,k,mesh_element(4,i))
! if (constitutive_phi1(l,j,i)==400*inRad) then
! call math_halton_ori()
! else
! call math_gauss()
! endif
! constitutive_phi1(l,j,i)=texture_Gauss(1,k,mesh_element(4,i))
! constitutive_phi(l,j,i)=texture_Gauss(2,k,mesh_element(4,i))
! constitutive_phi2(l,j,i)=texture_Gauss(3,k,mesh_element(4,i))
enddo
matID=mesh_element(3,i)
texID=mesh_element(4,i)
if (texture_ODFfile(texID)=='') then
multiplicity=texture_Ngrains(texID)/(texture_NGauss(texID)+texture_NFiber(texID)+texture_NRandom(texID))
do j=1,FE_Nips(mesh_element(2,i))
!* Gauss component
do k=1,multiplicity*texture_NGauss(texID),multiplicity
do l=k,k+multiplicity
constitutive_matID(l,j,i)=matID
constitutive_texID(l,j,i)=texID
constitutive_MatVolFrac(l,j,i)=1.0_pReal
constitutive_TexVolFrac(l,j,i)=texture_Gauss(6,k,texID)/multiplicity
constitutive_phi1(l,j,i)=texture_Gauss(1,k,texID)
constitutive_phi(l,j,i)=texture_Gauss(2,k,texID)
constitutive_phi2(l,j,i)=texture_Gauss(3,k,texID)
! if (constitutive_phi1(l,j,i)==400*inRad) then
! call math_halton_ori()
! else
! call math_gauss()
! endif
! constitutive_phi1(l,j,i)=texture_Gauss(1,k,mesh_element(4,i))
! constitutive_phi(l,j,i)=texture_Gauss(2,k,mesh_element(4,i))
! constitutive_phi2(l,j,i)=texture_Gauss(3,k,mesh_element(4,i))
enddo
enddo
!* Fiber component
do k=1+texture_NGauss(texID),multiplicity*texture_NFiber(texID)+texture_NGauss(texID),multiplicity
do l=k,k+multiplicity
constitutive_matID(l,j,i)=matID
constitutive_texID(l,j,i)=texID
constitutive_MatVolFrac(l,j,i)=1.0_pReal
constitutive_TexVolFrac(l,j,i)=texture_Fiber(6,k,texID)/multiplicity
! constitutive_phi1(l,j,i)=texture_Fiber(1,k,mesh_element(4,i))
! constitutive_phi(l,j,i)=texture_Fiber(2,k,mesh_element(4,i))
! constitutive_phi2(l,j,i)=texture_Fiber(3,k,mesh_element(4,i))
enddo
enddo
!* Random component
do k=1+texture_NGauss(texID)+texture_NFiber(texID),constitutive_Ngrains(j,i),multiplicity
do l=k,k+multiplicity
constitutive_matID(l,j,i)=matID
constitutive_texID(l,j,i)=texID
constitutive_MatVolFrac(l,j,i)=1.0_pReal
constitutive_TexVolFrac(l,j,i)=(1.0_pReal-texture_Gauss(5,k,texID)-texture_Fiber(6,k,texID))/multiplicity
! constitutive_phi1(l,j,i)=texture_Fiber(1,k,mesh_element(4,i))
! constitutive_phi(l,j,i)=texture_Fiber(2,k,mesh_element(4,i))
! constitutive_phi2(l,j,i)=texture_Fiber(3,k,mesh_element(4,i))
enddo
enddo
enddo
!* Fiber component
do k=1+texture_NGauss(mesh_element(4,i)),constitutive_Ngrains(i,j),mul
do l=k,k+mul
constitutive_matID(l,j,i)=mesh_element(3,i)
constitutive_texID(l,j,i)=mesh_element(4,i)
constitutive_MatVolFrac(l,j,i)=1.0_pReal
constitutive_TexVolFrac(l,j,i)=texture_Fiber(6,k,mesh_element(4,i))/mul
! constitutive_phi1(l,j,i)=texture_Fiber(1,k,mesh_element(4,i))
! constitutive_phi(l,j,i)=texture_Fiber(2,k,mesh_element(4,i))
! constitutive_phi2(l,j,i)=texture_Fiber(3,k,mesh_element(4,i))
enddo
enddo
enddo
endif
enddo
!* Initialization of state variables
!do l=1,material_Nstatevars(k,j,i)
! constitutive_state_old(l,k,j,i)=material_s0_slip(constitutive_matID(k,j,i))
@ -914,21 +948,17 @@ Tstar_v(6)=Tstar_v_m(6)/dsqrt(2.0_pReal)
Lp=0.0_pReal
do i=1,material_Nslip(matID)
tau_slip(i)=dot_product(Tstar_v,constitutive_Sslip_v(:,i,material_CrystalStructure(matID)))
gdot_slip(i)=material_gdot0_slip(matID)*(abs(tau_slip(i))/constitutive_state_new(i,ipc,ip,el))**material_n_slip(matID)* &
sign(1.0_pReal,tau_slip(i))
gdot_slip(i)=material_gdot0_slip(matID)*(abs(tau_slip(i))/constitutive_state_new(i,ipc,ip,el))**material_n_slip(matID)*sign(1.0_pReal,tau_slip(i))
Lp=Lp+gdot_slip(i)*constitutive_Sslip(:,:,i,material_CrystalStructure(matID))
enddo
!* Calculation of the tangent of Lp
dLp_dTstar=0.0_pReal
do i=1,material_Nslip(matID)
dgdot_dtauslip(i)=material_gdot0_slip(matID)* &
(abs(tau_slip(i))/constitutive_state_new(i,ipc,ip,el))**(material_n_slip(matID)-1.0_pReal)* &
material_n_slip(matID)/constitutive_state_new(i,ipc,ip,el)
dgdot_dtauslip(i)=material_gdot0_slip(matID)*(abs(tau_slip(i))/constitutive_state_new(i,ipc,ip,el))**(material_n_slip(matID)-1.0_pReal)*material_n_slip(matID)/constitutive_state_new(i,ipc,ip,el)
forall (k=1:3,l=1:3,m=1:3,n=1:3)
dLp_dTstar(k,l,m,n)=dLp_dTstar(k,l,m,n)+constitutive_Sslip(k,l,i,material_CrystalStructure(matID))* &
constitutive_Sslip(m,n,i,material_CrystalStructure(matID))*dgdot_dtauslip(i)
end forall
dLp_dTstar(k,l,m,n)=dLp_dTstar(k,l,m,n)+constitutive_Sslip(k,l,i,material_CrystalStructure(matID))*constitutive_Sslip(m,n,i,material_CrystalStructure(matID))*dgdot_dtauslip(i)
endforall
enddo
return
@ -970,17 +1000,12 @@ Tstar_v(6)=Tstar_v_m(6)/dsqrt(2.0_pReal)
!* Self-Hardening of each system
do i=1,constitutive_Nstatevars(ipc,ip,el)
tau_slip(i)=dot_product(Tstar_v,constitutive_Sslip_v(:,i,material_CrystalStructure(matID)))
gdot_slip(i)=material_gdot0_slip(matID)*(abs(tau_slip(i))/constitutive_state_new(i,ipc,ip,el))**material_n_slip(matID)* &
sign(1.0_pReal,tau_slip(i))
self_hardening(i)=material_h0(matID)*(1.0_pReal-constitutive_state_new(i,ipc,ip,el)/material_s_sat(matID))**material_w0(matID)* &
abs(gdot_slip(i))
gdot_slip(i)=material_gdot0_slip(matID)*(abs(tau_slip(i))/constitutive_state_new(i,ipc,ip,el))**material_n_slip(matID)*sign(1.0_pReal,tau_slip(i))
self_hardening(i)=material_h0(matID)*(1.0_pReal-constitutive_state_new(i,ipc,ip,el)/material_s_sat(matID))**material_w0(matID)*abs(gdot_slip(i))
enddo
!* Hardening for all systems
constitutive_DotState = matmul(constitutive_HardeningMatrix(1:material_Nslip(matID), &
1:material_Nslip(matID), &
material_CrystalStructure(matID)), &
self_hardening)
constitutive_DotState=matmul(constitutive_HardeningMatrix(1:material_Nslip(matID),1:material_Nslip(matID),material_CrystalStructure(matID)),self_hardening)
return
end function