DAMASK_EICMD/trunk/constitutive.f90

707 lines
29 KiB
Fortran
Raw Normal View History

!************************************
!* Module: CONSTITUTIVE *
!************************************
!* contains: *
!* - constitutive equations *
!* - Schmid matrices calculation *
!* - Hardening matrices definition *
!* - Parameters definition *
!* - orientations? *
!************************************
MODULE constitutive
!*** Include other modules ***
use prec, only: pReal,pInt
implicit none
! QUESTION: would it be wise to outsource these to _constitutive_ ?? YES!
! *** Slip resistances at (t=t0) and (t=t1) ***
! allocate(constitutive_state_old(constitutive_Nstatevars,constitutive_maxNgrains,mesh_maxNips,mesh_NcpElems))
!allocate(constitutive_state_new(constitutive_Nstatevars,constitutive_maxNgrains,mesh_maxNips,mesh_NcpElems))
!state_tauc_slip_old = 0.0_pReal
!state_tauc_slip_new = 0.0_pReal
! *** Transformation to get the MARC order ***
! *** 11,22,33,12,23,13 ***
! MISSING this should be outsourced to FEM-spec
!temp=Cslip_66(4,:)
!Cslip_66(4,:)=Cslip_66(6,:)
!Cslip_66(6,:)=Cslip_66(5,:)
!Cslip_66(5,:)=temp
!temp=Cslip_66(:,4)
!Cslip_66(:,4)=2.0d0*Cslip_66(:,6)
!Cslip_66(:,6)=2.0d0*Cslip_66(:,5)
!Cslip_66(:,5)=2.0d0*temp
!***********************************************
!* Definition of crystal structures properties *
!***********************************************
!* Number of crystal structures (1-FCC,2-BCC,3-HCP)
integer(pInt), parameter :: constitutive_MaxCrystalStructure = 3
!* Total number of slip systems per crystal structure
!* (as to be changed according the definition of slip systems)
integer(pInt), dimension(constitutive_MaxCrystalStructure), parameter :: constitutive_MaxNslipOfStructure = &
reshape((/12,48,12/),(/constitutive_MaxCrystalStructure/))
!* Maximum number of slip systems over crystal structures
integer(pInt), parameter :: constitutive_MaxMaxNslipOfStructure = 48
!* Slip direction, slip normales and Schmid matrices
real(pReal), dimension(3,3,constitutive_MaxMaxNslipOfStructure,constitutive_MaxCrystalStructure) :: constitutive_Sslip
real(pReal), dimension(6,constitutive_MaxMaxNslipOfStructure,constitutive_MaxCrystalStructure) :: constitutive_Sslip_v
real(pReal), dimension(3,constitutive_MaxMaxNslipOfStructure,constitutive_MaxCrystalStructure) :: constitutive_sn
real(pReal), dimension(3,constitutive_MaxMaxNslipOfStructure,constitutive_MaxCrystalStructure) :: constitutive_sd
!*** Slip systems for FCC structures (1) ***
!* System {111}<110> Sort according Eisenlohr&Hantcherli
data constitutive_sd(:, 1,1)/ 0, 1,-1/ ; data constitutive_sn(:, 1,1)/ 1, 1, 1/
data constitutive_sd(:, 2,1)/-1, 0, 1/ ; data constitutive_sn(:, 2,1)/ 1, 1, 1/
data constitutive_sd(:, 3,1)/ 1,-1, 0/ ; data constitutive_sn(:, 3,1)/ 1, 1, 1/
data constitutive_sd(:, 4,1)/ 0,-1,-1/ ; data constitutive_sn(:, 4,1)/-1,-1, 1/
data constitutive_sd(:, 5,1)/ 1, 0, 1/ ; data constitutive_sn(:, 5,1)/-1,-1, 1/
data constitutive_sd(:, 6,1)/-1, 1, 0/ ; data constitutive_sn(:, 6,1)/-1,-1, 1/
data constitutive_sd(:, 7,1)/ 0,-1, 1/ ; data constitutive_sn(:, 7,1)/ 1,-1,-1/
data constitutive_sd(:, 8,1)/-1, 0,-1/ ; data constitutive_sn(:, 8,1)/ 1,-1,-1/
data constitutive_sd(:, 9,1)/ 1, 1, 0/ ; data constitutive_sn(:, 9,1)/ 1,-1,-1/
data constitutive_sd(:,10,1)/ 0, 1, 1/ ; data constitutive_sn(:,10,1)/-1, 1,-1/
data constitutive_sd(:,11,1)/ 1, 0,-1/ ; data constitutive_sn(:,11,1)/-1, 1,-1/
data constitutive_sd(:,12,1)/-1,-1, 0/ ; data constitutive_sn(:,12,1)/-1, 1,-1/
!*** Slip systems for BCC structures (2) ***
!* System {110}<111>
!* Sort?
data constitutive_sd(:, 1,2)/ 1,-1, 1/ ; data constitutive_sn(:, 1,2)/ 0, 1, 1/
data constitutive_sd(:, 2,2)/-1,-1, 1/ ; data constitutive_sn(:, 2,2)/ 0, 1, 1/
data constitutive_sd(:, 3,2)/ 1, 1, 1/ ; data constitutive_sn(:, 3,2)/ 0,-1, 1/
data constitutive_sd(:, 4,2)/-1, 1, 1/ ; data constitutive_sn(:, 4,2)/ 0,-1, 1/
data constitutive_sd(:, 5,2)/-1, 1, 1/ ; data constitutive_sn(:, 5,2)/ 1, 0, 1/
data constitutive_sd(:, 6,2)/-1,-1, 1/ ; data constitutive_sn(:, 6,2)/ 1, 0, 1/
data constitutive_sd(:, 7,2)/ 1, 1, 1/ ; data constitutive_sn(:, 7,2)/-1, 0, 1/
data constitutive_sd(:, 8,2)/ 1,-1, 1/ ; data constitutive_sn(:, 8,2)/-1, 0, 1/
data constitutive_sd(:, 9,2)/-1, 1, 1/ ; data constitutive_sn(:, 9,2)/ 1, 1, 0/
data constitutive_sd(:,10,2)/-1, 1,-1/ ; data constitutive_sn(:,10,2)/ 1, 1, 0/
data constitutive_sd(:,11,2)/ 1, 1, 1/ ; data constitutive_sn(:,11,2)/-1, 1, 0/
data constitutive_sd(:,12,2)/ 1, 1,-1/ ; data constitutive_sn(:,12,2)/-1, 1, 0/
!* System {112}<111>
!* Sort?
data constitutive_sd(:,13,2)/-1, 1, 1/ ; data constitutive_sn(:,13,2)/ 2, 1, 1/
data constitutive_sd(:,14,2)/ 1, 1, 1/ ; data constitutive_sn(:,14,2)/-2, 1, 1/
data constitutive_sd(:,15,2)/ 1, 1,-1/ ; data constitutive_sn(:,15,2)/ 2,-1, 1/
data constitutive_sd(:,16,2)/ 1,-1, 1/ ; data constitutive_sn(:,16,2)/ 2, 1,-1/
data constitutive_sd(:,17,2)/ 1,-1, 1/ ; data constitutive_sn(:,17,2)/ 1, 2, 1/
data constitutive_sd(:,18,2)/ 1, 1,-1/ ; data constitutive_sn(:,18,2)/-1, 2, 1/
data constitutive_sd(:,19,2)/ 1, 1, 1/ ; data constitutive_sn(:,19,2)/ 1,-2, 1/
data constitutive_sd(:,20,2)/-1, 1, 1/ ; data constitutive_sn(:,20,2)/ 1, 2,-1/
data constitutive_sd(:,21,2)/ 1, 1,-1/ ; data constitutive_sn(:,21,2)/ 1, 1, 2/
data constitutive_sd(:,22,2)/ 1,-1, 1/ ; data constitutive_sn(:,22,2)/-1, 1, 2/
data constitutive_sd(:,23,2)/-1, 1, 1/ ; data constitutive_sn(:,23,2)/ 1,-1, 2/
data constitutive_sd(:,24,2)/ 1, 1, 1/ ; data constitutive_sn(:,24,2)/ 1, 1,-2/
!* System {123}<111>
!* Sort?
data constitutive_sd(:,25,2)/ 1, 1,-1/ ; data constitutive_sn(:,25,2)/ 1, 2, 3/
data constitutive_sd(:,26,2)/ 1,-1, 1/ ; data constitutive_sn(:,26,2)/-1, 2, 3/
data constitutive_sd(:,27,2)/-1, 1, 1/ ; data constitutive_sn(:,27,2)/ 1,-2, 3/
data constitutive_sd(:,28,2)/ 1, 1, 1/ ; data constitutive_sn(:,28,2)/ 1, 2,-3/
data constitutive_sd(:,29,2)/ 1,-1, 1/ ; data constitutive_sn(:,29,2)/ 1, 3, 2/
data constitutive_sd(:,30,2)/ 1, 1,-1/ ; data constitutive_sn(:,30,2)/-1, 3, 2/
data constitutive_sd(:,31,2)/ 1, 1, 1/ ; data constitutive_sn(:,31,2)/ 1,-3, 2/
data constitutive_sd(:,32,2)/-1, 1, 1/ ; data constitutive_sn(:,32,2)/ 1, 3,-2/
data constitutive_sd(:,33,2)/ 1, 1,-1/ ; data constitutive_sn(:,33,2)/ 2, 1, 3/
data constitutive_sd(:,34,2)/ 1,-1, 1/ ; data constitutive_sn(:,34,2)/-2, 1, 3/
data constitutive_sd(:,35,2)/-1, 1, 1/ ; data constitutive_sn(:,35,2)/ 2,-1, 3/
data constitutive_sd(:,36,2)/ 1, 1, 1/ ; data constitutive_sn(:,36,2)/ 2, 1,-3/
data constitutive_sd(:,37,2)/ 1,-1, 1/ ; data constitutive_sn(:,37,2)/ 2, 3, 1/
data constitutive_sd(:,38,2)/ 1, 1,-1/ ; data constitutive_sn(:,38,2)/-2, 3, 1/
data constitutive_sd(:,39,2)/ 1, 1, 1/ ; data constitutive_sn(:,39,2)/ 2,-3, 1/
data constitutive_sd(:,40,2)/-1, 1, 1/ ; data constitutive_sn(:,40,2)/ 2, 3,-1/
data constitutive_sd(:,41,2)/-1, 1, 1/ ; data constitutive_sn(:,41,2)/ 3, 1, 2/
data constitutive_sd(:,42,2)/ 1, 1, 1/ ; data constitutive_sn(:,42,2)/-3, 1, 2/
data constitutive_sd(:,43,2)/ 1, 1,-1/ ; data constitutive_sn(:,43,2)/ 3,-1, 2/
data constitutive_sd(:,44,2)/ 1,-1, 1/ ; data constitutive_sn(:,44,2)/ 3, 1,-2/
data constitutive_sd(:,45,2)/-1, 1, 1/ ; data constitutive_sn(:,45,2)/ 3, 2, 1/
data constitutive_sd(:,46,2)/ 1, 1, 1/ ; data constitutive_sn(:,46,2)/-3, 2, 1/
data constitutive_sd(:,47,2)/ 1, 1,-1/ ; data constitutive_sn(:,47,2)/ 3,-2, 1/
data constitutive_sd(:,48,2)/ 1,-1, 1/ ; data constitutive_sn(:,48,2)/ 3, 2,-1/
!*** Slip systems for HCP structures (3) ***
!* Basal systems {0001}<1120> (independent of c/a-ratio)
!* 1- (0 0 0 1)[-2 1 1 0]
!* 2- (0 0 0 1)[ 1 -2 1 0]
!* 3- (0 0 0 1)[ 1 1 -2 0]
!* Plane (hkil)->(hkl)
!* Direction [uvtw]->[(u-t) (v-t) w]
!* Automatical transformation from Bravais to Miller
!* not done for the moment
!* Sort?
data constitutive_sd(:, 1,3)/-1, 0, 0/ ; data constitutive_sn(:, 1,3)/ 0, 0, 1/
data constitutive_sd(:, 2,3)/ 0,-1, 0/ ; data constitutive_sn(:, 2,3)/ 0, 0, 1/
data constitutive_sd(:, 3,3)/ 1, 1, 0/ ; data constitutive_sn(:, 3,3)/ 0, 0, 1/
!* 1st type prismatic systems {1010}<1120> (independent of c/a-ratio)
!* 1- ( 0 1 -1 0)[-2 1 1 0]
!* 2- ( 1 0 -1 0)[ 1 -2 1 0]
!* 3- (-1 1 0 0)[ 1 1 -2 0]
!* Sort?
data constitutive_sd(:, 4,3)/-1, 0, 0/ ; data constitutive_sn(:, 4,3)/ 0, 1, 0/
data constitutive_sd(:, 5,3)/ 0,-1, 0/ ; data constitutive_sn(:, 5,3)/ 1, 0, 0/
data constitutive_sd(:, 6,3)/ 1, 1, 0/ ; data constitutive_sn(:, 6,3)/-1, 1, 0/
!* 1st type 1st order pyramidal systems {1011}<1120>
!* plane normales depend on the c/a-ratio
!* 1- ( 0 -1 1 1)[-2 1 1 0]
!* 2- ( 0 1 -1 1)[-2 1 1 0]
!* 3- (-1 0 1 1)[ 1 -2 1 0]
!* 4- ( 1 0 -1 1)[ 1 -2 1 0]
!* 5- (-1 1 0 1)[ 1 1 -2 0]
!* 6- ( 1 -1 0 1)[ 1 1 -2 0]
!* Sort?
data constitutive_sd(:, 7,3)/-1, 0, 0/ ; data constitutive_sn(:, 7,3)/ 0,-1, 1/
data constitutive_sd(:, 8,3)/ 0,-1, 0/ ; data constitutive_sn(:, 8,3)/ 0, 1, 1/
data constitutive_sd(:, 9,3)/ 1, 1, 0/ ; data constitutive_sn(:, 9,3)/-1, 0, 1/
data constitutive_sd(:,10,3)/-1, 0, 0/ ; data constitutive_sn(:,10,3)/ 1, 0, 1/
data constitutive_sd(:,11,3)/ 0,-1, 0/ ; data constitutive_sn(:,11,3)/-1, 1, 1/
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_hardening_matrix
real(pReal), parameter :: constitutive_latent_hardening=1.4
!*************************************
!* Definition of material properties *
!*************************************
!* Number of materials
integer(pInt) materials_maxN
!* Crystal structure and number of selected slip systems per material
integer(pInt), dimension(:) , allocatable :: materials_CrystalStructure
integer(pInt), dimension(:) , allocatable :: materials_Nslip
!* Maximum number of selected slip systems over materials
integer(pInt), allocatable :: materials_MaxNslip
!* Elastic constants and matrices
real(pReal), dimension(:) , allocatable :: materials_C11
real(pReal), dimension(:) , allocatable :: materials_C12
real(pReal), dimension(:) , allocatable :: materials_C13
real(pReal), dimension(:) , allocatable :: materials_C33
real(pReal), dimension(:) , allocatable :: materials_C44
real(pReal), dimension(:,:,:), allocatable :: materials_Cslip_66
! NB: Cslip_66(1:6,1:6,number of materials)
!* Visco-plastic material parameters
real(pReal), dimension(:) , allocatable :: materials_s0_slip
real(pReal), dimension(:) , allocatable :: materials_gdot0_slip
real(pReal), dimension(:) , allocatable :: materials_n_slip
real(pReal), dimension(:) , allocatable :: materials_h0
real(pReal), dimension(:) , allocatable :: materials_s_sat
real(pReal), dimension(:) , allocatable :: materials_w0
! NB: Parameters(number of materials)
!************************************
!* Definition of texture properties *
!************************************
!* Number of textures
integer(pInt) textures_maxN
!* Textures definition
character(len=80), dimension(:), allocatable :: textures_ODFfile
character(len=80), dimension(:), allocatable :: textures_symmetry
integer(pInt), dimension(:) , allocatable :: textures_Ngrains
! NB: symmetry(number of texture)
CONTAINS
!****************************************
!* - constitutive_init
!* - constitutive_SchmidMatrices
!* - constitutive_HardeningMatrices
!* - constitutive_CountSections
!* - constitutive_Parse_UnknownPart
!* - constitutive_Parse_MaterialPart
!* - constitutive_Parse_TexturePart
!* - constitutive_Parse_MatTexDat
!* - constitutive_Lp
!* - constitutive_TangentLp
!* - consistutive_DotState
!****************************************
subroutine constitutive_init()
!**************************************
!* Module initialization *
!**************************************
call constitutive_calc_SchmidMatrices()
call constitutive_calc_HardeningMatrices()
call constitutive_parse_MatTexDat('materials_textures.mpie')
end subroutine
subroutine constitutive_SchmidMatrices()
!**************************************
!* Calculation of Schmid matrices *
!**************************************
use prec, only: pReal,pInt
implicit none
!* Definition of variables
integer(pInt) i,j,k,l
real(pReal) invNorm
!* Iteration over the crystal structures
do l=1,3
!* Iteration over the systems
do k=1,constitutive_MaxNslipOfStructure(l)
!* Defintion of Schmid matrix
forall (i=1:3,j=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)))
constitutive_Sslip(:,:,k,l)=constitutive_Sslip(:,:,k,l)*invNorm
!* Vectorization of normalized Schmid matrix
!* according MARC component order 11,22,33,12,23,13
constitutive_Sslip_v(1,k,l)=constitutive_Sslip(1,1,k,l)
constitutive_Sslip_v(2,k,l)=constitutive_Sslip(2,2,k,l)
constitutive_Sslip_v(3,k,l)=constitutive_Sslip(3,3,k,l)
constitutive_Sslip_v(4,k,l)=constitutive_Sslip(1,2,k,l)+constitutive_Sslip(2,1,k,l)
constitutive_Sslip_v(5,k,l)=constitutive_Sslip(2,3,k,l)+constitutive_Sslip(3,3,k,l)
constitutive_Sslip_v(6,k,l)=constitutive_Sslip(1,3,k,l)+constitutive_Sslip(3,1,k,l)
enddo
enddo
end subroutine
subroutine constitutive_HardeningMatrices()
!****************************************
!* Hardening matrix (see Kalidindi) *
!****************************************
use prec, only: pReal,pInt
implicit none
!* Definition of variables
integer(pInt) i,j,k,l
!* Initialization of the hardening matrix
constitutive_hardening_matrix=constitutive_latent_hardening
!* Iteration over the crystal structures
do l=1,3
select case(l)
!* Hardening matrix for FCC structures
case (1)
do k=1,10,3
forall (i=1:3,j=1:3)
constitutive_hardening_matrix(k-1+i,k-1+j,l)=1.0_pReal
endforall
enddo
!* Hardening matrix for BCC structures
case (2)
do k=1,11,2
forall (i=1:2,j=1:2)
constitutive_hardening_matrix(k-1+i,k-1+j,l)=1.0_pReal
endforall
enddo
do k=13,48
constitutive_hardening_matrix(k,k,l)=1.0_pReal
enddo
!* Hardening matrix for HCP structures
case (3)
forall (i=1:3,j=1:3)
constitutive_hardening_matrix(i,j,l)=1.0_pReal
endforall
do k=4,12
constitutive_hardening_matrix(k,k,l)=1.0_pReal
enddo
end select
enddo
end subroutine
subroutine constitutive_CountSections(file,count,part)
!*********************************************************************
!* This subroutine reads a "part" from the input file until the next *
!* part is reached and counts the number of "sections" in the part *
!* INPUT: *
!* - file : file ID *
!* OUTPUT: *
!* - part : name of the next "part" *
!* - count : number of sections inside the current "part" *
!*********************************************************************
use prec, only: pInt
use IO, only: IO_stringPos,IO_stringValue,IO_lc
implicit none
!* Definition of variables
character(len=80) part,line,tag
integer(pInt) file,count,pos
integer(pInt), dimension(3) :: positions
count=0
part=''
do while(.true.)
read(unit=file,fmt='(a80)',END=100) line
positions=IO_stringPos(line,1)
tag=IO_lc(IO_stringValue(line,positions,1))
if (tag(1:1)=='<' .and. tag(len(tag):len(tag)=='>') then
part=tag(2:len(tag)-1)
exit
elseif (tag(1:1)=='[' .and. tag(len(tag):len(tag)==']') then
count=count+1
endif
enddo
100 return
end subroutine
character(len=80) function constitutive_Parse_UnknownPart(file)
!*********************************************************************
!* This function reads a unknown "part" from the input file until *
!* the next part is reached *
!* INPUT: *
!* - file : file ID *
!*********************************************************************
use prec, only: pInt
use IO, only: IO_stringPos,IO_stringValue,IO_lc
implicit none
!* Definition of variables
character(len=80) line,tag
integer(pInt), parameter :: maxNchunks = 1
integer(pInt) file
integer(pInt), dimension(1+2*maxNchunks) :: positions
constitutive_parse_unknownPart = ''
do while(.true.)
read(unit=file,fmt='(a80)',end=100) line
positions = IO_stringPos(line,maxNchunks)
tag = IO_lc(IO_stringValue(line,positions,1))
if (tag(1:1)=='<' .and. tag(len(tag):len(tag)=='>') then
constitutive_parse_unknownPart = tag(2:len(tag)-1)
exit
endif
enddo
100 return
end function
character(len=80) function constitutive_Parse_MaterialPart(file)
!*********************************************************************
!* This function reads a material "part" from the input file until *
!* the next part is reached *
!* INPUT: *
!* - file : file ID *
!*********************************************************************
use prec, only: pInt
use IO, only: IO_stringPos,IO_stringValue,IO_lc
implicit none
!* Definition of variables
character(len=80) line,tag
integer(pInt), parameter :: maxNchunks = 2! may be more than 2 chunks ..?
integer(pInt) file,section
integer(pInt), dimension(1+2*maxNchunks) :: positions
section = 0
constitutive_parse_materialPart = ''
do while(.true.)
read(unit=file,fmt='(a80)',end=100) line
positions = IO_stringPos(line,maxNchunks) ! parse leading chunks
tag = IO_lc(IO_stringValue(line,positions,1))
if (tag(1:1)=='#') then ! skip comment line
cycle
elseif (tag(1:1)=='<' .and. tag(len(tag):len(tag)=='>') then
constitutive_parse_materialPart = tag(2:len(tag)-1)
exit
elseif (tag(1:1)=='[' .and. tag(len(tag):len(tag)==']') then
section = section+1
else
if (section>0) then
select case(tag)
case ('crystal_structure')
materials_CrystalStructure(section)=IO_intValue(line,positions,2)
case ('nslip')
materials_Nslip(section)=IO_intValue(line,positions,2)
case ('C11')
materials_C11(section)=IO_floatValue(line,positions,2)
case ('C12')
materials_C12(section)=IO_floatValue(line,positions,2)
case ('C13')
materials_C13(section)=IO_floatValue(line,positions,2)
case ('C33')
materials_C33(section)=IO_floatValue(line,positions,2)
case ('C44')
materials_C44(section)=IO_floatValue(line,positions,2)
case ('s0_slip')
materials_s0_slip(section)=IO_floatValue(line,positions,2)
case ('gdot0_slip')
materials_gdot0_slip(section)=IO_floatValue(line,positions,2)
case ('n_slip')
materials_n_slip(section)=IO_floatValue(line,positions,2)
case ('h0')
materials_h0(section)=IO_floatValue(line,positions,2)
case ('s_sat')
materials_s_sat(section)=IO_floatValue(line,positions,2)
case ('w0')
materials_w0(section)=IO_floatValue(line,positions,2)
case default
write(6,*) 'Unknown material parameter ',line
end select
endif
endif
enddo
100 return
end function
character(len=80) function constitutive_Parse_TexturePart(file)
!*********************************************************************
!* This function reads a texture "part" from the input file until *
!* the next part is reached *
!* INPUT: *
!* - file : file ID *
!*********************************************************************
use prec, only: pInt
use IO, only: IO_stringPos,IO_stringValue,IO_lc
implicit none
!* Definition of variables
character(len=80) line,tag
integer(pInt), parameter :: maxNchunks = 10 ! may be more than 10 chunks ..?
integer(pInt) file,pos,section
integer(pInt), dimension(1+2*maxNchunks) :: positions
section = 0
constitutive_parse_texturePart = ''
do while(.true.)
read(unit=file,fmt='(a80)',end=100) line
positions = IO_stringPos(line,maxNchunks) ! parse leading chunks
tag = IO_lc(IO_stringValue(line,positions,1))
if (tag(1:1)=='#') then ! skip comment line
cycle
elseif (tag(1:1)=='<' .and. tag(len(tag):len(tag)=='>') then
constitutive_parse_texturePart = tag(2:len(tag)-1)
exit
elseif (tag(1:1)=='[' .and. tag(len(tag):len(tag)==']') then
section = section+1
else
if (section>0) then
select case(tag)
case ('hybridIA')
textures_ODFfile(section)=IO_stringValue(line,positions,2)
case ('gauss')
case ('fiber')
case ('ngrains')
textures_Ngrains(section)=IO_intValue(line,positions,2)
case ('symmetry')
textures_symmetry(section)=IO_stringValue(line,positions,2)
case default
write(6,*) 'Unknown texture parameter ',line
end select
endif
endif
enddo
100 return
end function
subroutine constitutive_Parse_MatTexDat(filename)
!*********************************************************************
!* This function reads the material and texture input file *
!* INPUT: *
!* - filename : name of input file *
!*********************************************************************
use prec, only: pReal,pInt
use IO
implicit none
!* Definition of variables
character(len=*) filename
character(len=80) part,formerPart
integer(pInt) sectionCount,i,j,m
!* Open input file
open(200,FILE=filename,ACTION='READ',STATUS='OLD',ERR=100)
!* First reading: number of materials and textures
!* Arrays allocation
part = '_dummy_'
do while (part/='')
formerPart = part
call constitutive_countSections(200,sectionCount,part)
select case (formerPart)
case ('materials')
materials_maxN = sectionCount
case ('textures')
textures_maxN = sectionCount
end select
enddo
allocate(textures_ODFfile(textures_maxN)) ; textures_ODFfile=''
allocate(textures_Ngrains(textures_maxN)) ; textures_Ngrains=0_pInt
allocate(textures_symmetry(textures_maxN)) ; textures_symmetry=''
allocate(materials_CrystalStructure(materials_maxN)) ; materials_CrystalStructure=0_pInt
allocate(materials_Nslip(materials_maxN)) ; materials_Nslip=0_pInt
allocate(materials_C11(materials_maxN)) ; materials_C11=0.0_pReal
allocate(materials_C12(materials_maxN)) ; materials_C12=0.0_pReal
allocate(materials_C13(materials_maxN)) ; materials_C13=0.0_pReal
allocate(materials_C33(materials_maxN)) ; materials_C33=0.0_pReal
allocate(materials_C44(materials_maxN)) ; materials_C44=0.0_pReal
allocate(materials_s0_slip(materials_maxN)) ; materials_s0_slip=0.0_pReal
allocate(materials_gdot0_slip(materials_maxN)) ; materials_gdot0_slip=0.0_pReal
allocate(materials_n_slip(materials_maxN)) ; materials_n_slip=0.0_pReal
allocate(materials_h0(materials_maxN)) ; materials_h0=0.0_pReal
allocate(materials_s_sat(materials_maxN)) ; materials_s_sat=0.0_pReal
allocate(materials_w0(materials_maxN)) ; materials_w0=0.0_pReal
!* Second reading: materials and textures are stored
part = '_dummy_'
do while (part/='')
select case (part)
case ('materials')
part = constitutive_parse_materialPart(200)
case ('textures')
part = constitutive_parse_texturePart(200)
case default
part = constitutive_parse_unknownPart(200)
end select
end do
!* Close input file
close(200)
do m=1,material_maxN
material_Cslip_66(:,:,m) = 0.0_pReal
select case (material_crystal_structure)
case (1:2) ! cubic structure
do i=1,3
do j=1,3
material_Cslip_66(i,j,m) = C12
enddo
material_Cslip_66(i,i,m) = C11
material_Cslip_66(i+3,i+3,m) = C44
enddo
case (3) ! hcp structure MISSING correct
do i=1,3
do j=1,3
material_Cslip_66(i,j,m) = C12
enddo
material_Cslip_66(i,i,m) = C11
material_Cslip_66(i+3,i+3,m) = C44
enddo
end select
material_Cslip_3333(:,:,:,:,m) = math_66to3333(Cslip_66(:,:,m))
end do
! MISSING some consistency checks may be..?
return
100 call IO_error(110) ! corrupt matarials_textures file
end subroutine
real(pReal) function constitutive_Lp(Tstar_v,ipc,ip,el)
!*********************************************************************
!* This subroutine contains the constitutive equation for *
!* calculating the velocity gradient *
!* INPUT: *
!* - Tstar_v : 2nd Piola Kirchhoff stress tensor *
!* - ipc : component-ID of current integration point *
!* - ip : current integration point *
!* - el : current element *
!*********************************************************************
use prec, only: pReal,pInt
implicit none
!* Definition of variables
!integer(pInt) matID,i
!real(pReal) tau_slip(constitutive_Nslip(matID))
!real(pReal) tauc_slip(constitutive_Nslip(matID))
!real(pReal) gdot_slip(constitutive_Nslip(matID))
!real(pReal) dgdot_dtaucslip(constitutive_Nslip(matID))
!* Iteration over the systems
!do i=1,constitutive_Nslip(matID)
! gdot_slip(i)=constitutive_gdot0_slip(matID)*(abs(tau_slip(i))/tauc_slip(i))**constitutive_n_slip(matID)*sign(1.0_pReal,tau_slip(i))
! dgdot_dtaucslip(i)=constitutive_gdot0_slip(matID)*(abs(tau_slip(i))/tauc_slip(i))**(constitutive_n_slip(matID)-1.0_pReal)*constitutive_n_slip(matID)/tauc_slip(i)
!enddo
return
end function
subroutine constitutive_TangentLp(Tstar_v,ipc,ip,el,dLp_dTstar,dLpT_dTstar)
!*********************************************************************
!* This subroutine contains the constitutive equation for *
!* calculating the velocity gradient *
!* INPUT: *
!* - Tstar_v : 2nd Piola Kirchhoff stress tensor *
!* - ipc : component-ID of current integration point *
!* - ip : current integration point *
!* - el : current element *
!* OUTPUT: *
!* - dLp_dTstar : derivative of Lp *
!* - dLpT_dTstar : derivative of tranposed Lp *
!*********************************************************************
use prec, only: pReal,pInt
implicit none
!* Definition of variables
!integer(pInt) matID,i,j
!real(pReal) tauc_slip(constitutive_Nslip(matID))
!real(pReal) gdot_slip(constitutive_Nslip(matID))
!real(pReal) dtauc_slip(constitutive_Nslip(matID))
!real(pReal) self_hardening(constitutive_Nslip(matID))
!* Self-Hardening of each system
!do i=1,constitutive_Nslip(matID)
! self_hardening(i)=constitutive_h0(matID)*(1.0_pReal-tauc_slip(i)/constitutive_s_sat(matID))**constitutive_w0(matID)*abs(gdot_slip(i))
!enddo
!* Hardening for all systems
!i=constitutive_Nslip(matID)
!j=constitutive_crystal_structure(matID)
!dtauc_slip=matmul(constitutive_hardening_matrix(1:i,1:i,j),self_hardening)
return
end subroutine
real(pReal) function constitutive_DotState(Tstar_v,ipc,ip,el)
!*********************************************************************
!* This subroutine contains the constitutive equation for *
!* calculating the velocity gradient *
!* INPUT: *
!* - Tstar_v : 2nd Piola Kirchhoff stress tensor *
!* - ipc : component-ID of current integration point *
!* - ip : current integration point *
!* - el : current element *
!*********************************************************************
use prec, only: pReal,pInt
implicit none
!* Definition of variables
!integer(pInt) matID,i
!real(pReal) dt,Lp(3,3)
!real(pReal) tau_slip(constitutive_Nslip(matID))
!real(pReal) tauc_slip_new(constitutive_Nslip(matID))
!real(pReal) gdot_slip(constitutive_Nslip(matID))
!* Calculation of Lp
!Lp=0.0_pReal
!do i=1,constitutive_Nslip(matID)
! gdot_slip(i)=constitutive_gdot0_slip(matID)*(abs(tau_slip(i))/tauc_slip(i))**constitutive_n_slip(matID)*sign(1.0_pReal,tau_slip(i))
! Lp=Lp+gdot_slip(i)*constitutive_Sslip(:,:,i,constitutive_crystal_structure(matID))
!enddo
return
end function
END MODULE