convergence of state in crystallite is now tested as follows:

(state < relevant state) or (residuum < relative tolerance * state)
since the relevant value for the state variables depend on their nature and can vary by large scales (e.g. volume fraction: 1e-10, dislocation density: 1e5) it is not possible to set a unique value. instead the constitutive law has to decide what is relevant. therefore, all constitutive laws now read in parameters from the material.config that determine the values for relevantState [@luc: in dislobased law relevant State is for the moment generally set to 1e-200, so no additional parameters necessary in material.config. if you also want this feature, we can still implement it, no big deal]

- added sanity checks in constitutive_nonlocal.f90

- corrected coordinate transformation for backstress calculation in constitutive_nonlocal.f90

- corrected equations for evolution of dipole dislocation densities (athermal annihilation and formation by glide)
This commit is contained in:
Christoph Kords 2009-09-18 15:37:14 +00:00
parent ec524a1dd2
commit b09b2b17f3
8 changed files with 239 additions and 83 deletions

View File

@ -33,7 +33,7 @@ subroutine IO_init
write(6,*) '$Id$'
write(6,*)
return
end subroutine
endsubroutine
!********************************************************************
! open existing file to given unit
@ -854,10 +854,28 @@ endfunction
msg = 'Non-positive diffusion prefactor'
case (225)
msg = 'No slip systems specified'
case (226)
msg = 'Non-positive prefactor for dislocation velocity'
case (227)
msg = 'Non-positive prefactor for mean free path'
case (228)
msg = 'Non-positive minimum stable dipole distance'
case (229)
msg = 'Hardening interaction coefficients below one'
case (230)
msg = 'Non-positive atomic volume'
case (231)
msg = 'Non-positive prefactor for self-diffusion coefficient'
case (232)
msg = 'Non-positive activation energy for dislocation climb'
case (233)
msg = 'Non-positive relevant dislocation density'
case (240)
msg = 'Non-positive Taylor factor'
case (241)
msg = 'Non-positive hardening exponent'
case (242)
msg = 'Non-positive relevant slip resistance'
case (260)
msg = 'Non-positive relevant strain'
case (261)

View File

@ -17,7 +17,8 @@ MODULE constitutive
constitutive_partionedState0, & ! pointer array to microstructure at start of homogenization inc
constitutive_subState0, & ! pointer array to microstructure at start of crystallite inc
constitutive_state, & ! pointer array to current microstructure (end of converged time step)
constitutive_dotState ! pointer array to evolution of current microstructure
constitutive_dotState, & ! pointer array to evolution of current microstructure
constitutive_relevantState ! relevant state values
integer(pInt), dimension(:,:,:), allocatable :: constitutive_sizeDotState, & ! size of dotState array
constitutive_sizeState, & ! size of state array per grain
constitutive_sizePostResults ! size of postResults array per grain
@ -110,58 +111,74 @@ subroutine constitutive_init()
allocate(constitutive_subState0(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems))
allocate(constitutive_state(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems))
allocate(constitutive_dotState(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems))
allocate(constitutive_relevantState(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems))
allocate(constitutive_sizeDotState(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems)) ; constitutive_sizeDotState = 0_pInt
allocate(constitutive_sizeState(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems)) ; constitutive_sizeState = 0_pInt
allocate(constitutive_sizePostResults(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems)); constitutive_sizePostResults = 0_pInt
do e = 1,mesh_NcpElems ! loop over elements
myNgrains = homogenization_Ngrains(mesh_element(3,e))
do i = 1,FE_Nips(mesh_element(2,e)) ! loop over IPs
do g = 1,myNgrains ! loop over grains
debugger = (e == 1 .and. i == 1 .and. g == 1)
myInstance = phase_constitutionInstance(material_phase(g,i,e))
select case(phase_constitution(material_phase(g,i,e)))
select case(phase_constitution(material_phase(g,i,e)))
case (constitutive_j2_label)
allocate(constitutive_state0(g,i,e)%p(constitutive_j2_sizeState(myInstance)))
allocate(constitutive_partionedState0(g,i,e)%p(constitutive_j2_sizeState(myInstance)))
allocate(constitutive_subState0(g,i,e)%p(constitutive_j2_sizeState(myInstance)))
allocate(constitutive_state(g,i,e)%p(constitutive_j2_sizeState(myInstance)))
allocate(constitutive_relevantState(g,i,e)%p(constitutive_j2_sizeState(myInstance)))
allocate(constitutive_dotState(g,i,e)%p(constitutive_j2_sizeDotState(myInstance)))
constitutive_state0(g,i,e)%p = constitutive_j2_stateInit(myInstance)
constitutive_relevantState(g,i,e)%p = constitutive_j2_relevantState(myInstance)
constitutive_sizeState(g,i,e) = constitutive_j2_sizeState(myInstance)
constitutive_sizeDotState(g,i,e) = constitutive_j2_sizeDotState(myInstance)
constitutive_sizePostResults(g,i,e) = constitutive_j2_sizePostResults(myInstance)
case (constitutive_phenopowerlaw_label)
allocate(constitutive_state0(g,i,e)%p(constitutive_phenopowerlaw_sizeState(myInstance)))
allocate(constitutive_partionedState0(g,i,e)%p(constitutive_phenopowerlaw_sizeState(myInstance)))
allocate(constitutive_subState0(g,i,e)%p(constitutive_phenopowerlaw_sizeState(myInstance)))
allocate(constitutive_state(g,i,e)%p(constitutive_phenopowerlaw_sizeState(myInstance)))
allocate(constitutive_relevantState(g,i,e)%p(constitutive_phenopowerlaw_sizeState(myInstance)))
allocate(constitutive_dotState(g,i,e)%p(constitutive_phenopowerlaw_sizeDotState(myInstance)))
constitutive_state0(g,i,e)%p = constitutive_phenopowerlaw_stateInit(myInstance)
constitutive_relevantState(g,i,e)%p = constitutive_phenopowerlaw_relevantState(myInstance)
constitutive_sizeState(g,i,e) = constitutive_phenopowerlaw_sizeState(myInstance)
constitutive_sizeDotState(g,i,e) = constitutive_phenopowerlaw_sizeDotState(myInstance)
constitutive_sizePostResults(g,i,e) = constitutive_phenopowerlaw_sizePostResults(myInstance)
case (constitutive_dislobased_label)
allocate(constitutive_state0(g,i,e)%p(constitutive_dislobased_sizeState(myInstance)))
allocate(constitutive_partionedState0(g,i,e)%p(constitutive_dislobased_sizeState(myInstance)))
allocate(constitutive_subState0(g,i,e)%p(constitutive_dislobased_sizeState(myInstance)))
allocate(constitutive_state(g,i,e)%p(constitutive_dislobased_sizeState(myInstance)))
allocate(constitutive_relevantState(g,i,e)%p(constitutive_dislobased_sizeState(myInstance)))
allocate(constitutive_dotState(g,i,e)%p(constitutive_dislobased_sizeDotState(myInstance)))
constitutive_state0(g,i,e)%p = constitutive_dislobased_stateInit(myInstance)
constitutive_relevantState(g,i,e)%p = constitutive_dislobased_relevantState(myInstance)
constitutive_sizeState(g,i,e) = constitutive_dislobased_sizeState(myInstance)
constitutive_sizeDotState(g,i,e) = constitutive_dislobased_sizeDotState(myInstance)
constitutive_sizePostResults(g,i,e) = constitutive_dislobased_sizePostResults(myInstance)
case (constitutive_nonlocal_label)
allocate(constitutive_state0(g,i,e)%p(constitutive_nonlocal_sizeState(myInstance)))
allocate(constitutive_partionedState0(g,i,e)%p(constitutive_nonlocal_sizeState(myInstance)))
allocate(constitutive_subState0(g,i,e)%p(constitutive_nonlocal_sizeState(myInstance)))
allocate(constitutive_state(g,i,e)%p(constitutive_nonlocal_sizeState(myInstance)))
allocate(constitutive_relevantState(g,i,e)%p(constitutive_nonlocal_sizeState(myInstance)))
allocate(constitutive_dotState(g,i,e)%p(constitutive_nonlocal_sizeDotState(myInstance)))
constitutive_state0(g,i,e)%p = constitutive_nonlocal_stateInit(myInstance)
constitutive_relevantState(g,i,e)%p = constitutive_nonlocal_relevantState(myInstance)
constitutive_sizeState(g,i,e) = constitutive_nonlocal_sizeState(myInstance)
constitutive_sizeDotState(g,i,e) = constitutive_nonlocal_sizeDotState(myInstance)
constitutive_sizePostResults(g,i,e) = constitutive_nonlocal_sizePostResults(myInstance)
case default
call IO_error(200,material_phase(g,i,e)) ! unknown constitution
end select
constitutive_partionedState0(g,i,e)%p = constitutive_state0(g,i,e)%p
enddo
@ -180,6 +197,7 @@ subroutine constitutive_init()
write(6,'(a32,x,7(i5,x))') 'constitutive_partionedState0: ', shape(constitutive_partionedState0)
write(6,'(a32,x,7(i5,x))') 'constitutive_subState0: ', shape(constitutive_subState0)
write(6,'(a32,x,7(i5,x))') 'constitutive_state: ', shape(constitutive_state)
write(6,'(a32,x,7(i5,x))') 'constitutive_relevantState: ', shape(constitutive_relevantState)
write(6,'(a32,x,7(i5,x))') 'constitutive_dotState: ', shape(constitutive_dotState)
write(6,'(a32,x,7(i5,x))') 'constitutive_sizeState: ', shape(constitutive_sizeState)
write(6,'(a32,x,7(i5,x))') 'constitutive_sizeDotState: ', shape(constitutive_sizeDotState)

View File

@ -543,6 +543,29 @@ function constitutive_dislobased_stateInit(myInstance)
end function
!*********************************************************************
!* relevant microstructural state *
!*********************************************************************
pure function constitutive_dislobased_relevantState(myInstance)
use prec, only: pReal, &
pInt
implicit none
!*** input variables
integer(pInt), intent(in) :: myInstance ! number specifying the current instance of the constitution
!*** output variables
real(pReal), dimension(constitutive_dislobased_sizeState(myInstance)) :: &
constitutive_dislobased_relevantState ! relevant state values for the current instance of this constitution
!*** local variables
constitutive_dislobased_relevantState = 1.0e-200_pReal
endfunction
function constitutive_dislobased_homogenizedC(state,ipc,ip,el)
!*********************************************************************
!* calculates homogenized elacticity matrix *

View File

@ -45,6 +45,7 @@ MODULE constitutive_j2
real(pReal), dimension(:), allocatable :: constitutive_j2_h0
real(pReal), dimension(:), allocatable :: constitutive_j2_tausat
real(pReal), dimension(:), allocatable :: constitutive_j2_w0
real(pReal), dimension(:), allocatable :: constitutive_j2_relevantResistance
CONTAINS
@ -82,23 +83,22 @@ subroutine constitutive_j2_init(file)
maxNinstance = count(phase_constitution == constitutive_j2_label)
if (maxNinstance == 0) return
allocate(constitutive_j2_sizeDotState(maxNinstance)) ; constitutive_j2_sizeDotState = 0_pInt
allocate(constitutive_j2_sizeState(maxNinstance)) ; constitutive_j2_sizeState = 0_pInt
allocate(constitutive_j2_sizePostResults(maxNinstance)); constitutive_j2_sizePostResults = 0_pInt
allocate(constitutive_j2_sizePostResult(maxval(phase_Noutput), &
maxNinstance)) ; constitutive_j2_sizePostResult = 0_pInt
allocate(constitutive_j2_output(maxval(phase_Noutput), &
maxNinstance)) ; constitutive_j2_output = ''
allocate(constitutive_j2_C11(maxNinstance)) ; constitutive_j2_C11 = 0.0_pReal
allocate(constitutive_j2_C12(maxNinstance)) ; constitutive_j2_C12 = 0.0_pReal
allocate(constitutive_j2_Cslip_66(6,6,maxNinstance)) ; constitutive_j2_Cslip_66 = 0.0_pReal
allocate(constitutive_j2_fTaylor(maxNinstance)) ; constitutive_j2_fTaylor = 0.0_pReal
allocate(constitutive_j2_tau0(maxNinstance)) ; constitutive_j2_tau0 = 0.0_pReal
allocate(constitutive_j2_gdot0(maxNinstance)) ; constitutive_j2_gdot0 = 0.0_pReal
allocate(constitutive_j2_n(maxNinstance)) ; constitutive_j2_n = 0.0_pReal
allocate(constitutive_j2_h0(maxNinstance)) ; constitutive_j2_h0 = 0.0_pReal
allocate(constitutive_j2_tausat(maxNinstance)) ; constitutive_j2_tausat = 0.0_pReal
allocate(constitutive_j2_w0(maxNinstance)) ; constitutive_j2_w0 = 0.0_pReal
allocate(constitutive_j2_sizeDotState(maxNinstance)) ; constitutive_j2_sizeDotState = 0_pInt
allocate(constitutive_j2_sizeState(maxNinstance)) ; constitutive_j2_sizeState = 0_pInt
allocate(constitutive_j2_sizePostResults(maxNinstance)); constitutive_j2_sizePostResults = 0_pInt
allocate(constitutive_j2_sizePostResult(maxval(phase_Noutput), maxNinstance)); constitutive_j2_sizePostResult = 0_pInt
allocate(constitutive_j2_output(maxval(phase_Noutput), maxNinstance)) ; constitutive_j2_output = ''
allocate(constitutive_j2_C11(maxNinstance)) ; constitutive_j2_C11 = 0.0_pReal
allocate(constitutive_j2_C12(maxNinstance)) ; constitutive_j2_C12 = 0.0_pReal
allocate(constitutive_j2_Cslip_66(6,6,maxNinstance)) ; constitutive_j2_Cslip_66 = 0.0_pReal
allocate(constitutive_j2_fTaylor(maxNinstance)) ; constitutive_j2_fTaylor = 0.0_pReal
allocate(constitutive_j2_tau0(maxNinstance)) ; constitutive_j2_tau0 = 0.0_pReal
allocate(constitutive_j2_gdot0(maxNinstance)) ; constitutive_j2_gdot0 = 0.0_pReal
allocate(constitutive_j2_n(maxNinstance)) ; constitutive_j2_n = 0.0_pReal
allocate(constitutive_j2_h0(maxNinstance)) ; constitutive_j2_h0 = 0.0_pReal
allocate(constitutive_j2_tausat(maxNinstance)) ; constitutive_j2_tausat = 0.0_pReal
allocate(constitutive_j2_w0(maxNinstance)) ; constitutive_j2_w0 = 0.0_pReal
allocate(constitutive_j2_relevantResistance(maxNinstance)) ; constitutive_j2_relevantResistance = 0.0_pReal
rewind(file)
line = ''
@ -142,17 +142,20 @@ subroutine constitutive_j2_init(file)
constitutive_j2_w0(i) = IO_floatValue(line,positions,2)
case ('taylorfactor')
constitutive_j2_fTaylor(i) = IO_floatValue(line,positions,2)
case ('relevantresistance')
constitutive_j2_relevantResistance(i) = IO_floatValue(line,positions,2)
end select
endif
enddo
100 do i = 1,maxNinstance ! sanity checks
if (constitutive_j2_tau0(i) < 0.0_pReal) call IO_error(210)
if (constitutive_j2_gdot0(i) <= 0.0_pReal) call IO_error(211)
if (constitutive_j2_n(i) <= 0.0_pReal) call IO_error(212)
if (constitutive_j2_tausat(i) <= 0.0_pReal) call IO_error(213)
if (constitutive_j2_w0(i) <= 0.0_pReal) call IO_error(241)
if (constitutive_j2_fTaylor(i) <= 0.0_pReal) call IO_error(240)
if (constitutive_j2_tau0(i) < 0.0_pReal) call IO_error(210)
if (constitutive_j2_gdot0(i) <= 0.0_pReal) call IO_error(211)
if (constitutive_j2_n(i) <= 0.0_pReal) call IO_error(212)
if (constitutive_j2_tausat(i) <= 0.0_pReal) call IO_error(213)
if (constitutive_j2_w0(i) <= 0.0_pReal) call IO_error(241)
if (constitutive_j2_fTaylor(i) <= 0.0_pReal) call IO_error(240)
if (constitutive_j2_relevantResistance(i) <= 0.0_pReal) call IO_error(242)
enddo
do i = 1,maxNinstance
@ -209,6 +212,29 @@ pure function constitutive_j2_stateInit(myInstance)
endfunction
!*********************************************************************
!* relevant microstructural state *
!*********************************************************************
pure function constitutive_j2_relevantState(myInstance)
use prec, only: pReal, &
pInt
implicit none
!*** input variables
integer(pInt), intent(in) :: myInstance ! number specifying the current instance of the constitution
!*** output variables
real(pReal), dimension(constitutive_j2_sizeState(myInstance)) :: &
constitutive_j2_relevantState ! relevant state values for the current instance of this constitution
!*** local variables
constitutive_j2_relevantState = constitutive_j2_relevantResistance(myInstance)
endfunction
function constitutive_j2_homogenizedC(state,ipc,ip,el)
!*********************************************************************
!* homogenized elacticity matrix *

View File

@ -17,18 +17,16 @@ implicit none
!* Definition of parameters
character (len=*), parameter :: constitutive_nonlocal_label = 'nonlocal'
character(len=16), dimension(9), parameter :: constitutive_nonlocal_stateList = (/ 'rhoEdgePos ', &
'rhoEdgeNeg ', &
'rhoScrewPos ', &
'rhoScrewNeg ', &
'rhoEdgeDip ', &
'rhoScrewDip ', &
'rhoForest ', &
'tauSlipThreshold', &
'Tdislocation_v ' /) ! list of microstructural state variables
character(len=16), dimension(6), parameter :: constitutive_nonlocal_stateListBasic = constitutive_nonlocal_stateList(1:6) ! list of "basic" microstructural state variables that are independent from other state variables
character(len=16), dimension(3), parameter :: constitutive_nonlocal_stateListDependent = constitutive_nonlocal_stateList(7:9) ! list of microstructural state variables that depend on other state variables
real(pReal), parameter :: kB = 1.38e-23_pReal ! Physical parameter, Boltzmann constant in J/Kelvin
character(len=16), dimension(6), parameter :: constitutive_nonlocal_listBasicStates = (/'rhoEdgePos ', &
'rhoEdgeNeg ', &
'rhoScrewPos ', &
'rhoScrewNeg ', &
'rhoEdgeDip ', &
'rhoScrewDip ' /) ! list of "basic" microstructural state variables that are independent from other state variables
character(len=16), dimension(3), parameter :: constitutive_nonlocal_listDependentStates = (/'rhoForest ', &
'tauSlipThreshold', &
'Tdislocation_v ' /) ! list of microstructural state variables that depend on other state variables
real(pReal), parameter :: kB = 1.38e-23_pReal ! Physical parameter, Boltzmann constant in J/Kelvin
!* Definition of global variables
@ -54,8 +52,9 @@ real(pReal), dimension(:), allocatable :: constitutive_nonlocal_
constitutive_nonlocal_Gmod, & ! shear modulus
constitutive_nonlocal_nu, & ! poisson's ratio
constitutive_nonlocal_atomicVolume, & ! atomic volume
constitutive_nonlocal_D0, & !
constitutive_nonlocal_Qsd
constitutive_nonlocal_D0, & ! prefactor for self-diffusion coefficient
constitutive_nonlocal_Qsd, & ! activation energy for dislocation climb
constitutive_nonlocal_relevantRho ! dislocation density considered relevant
real(pReal), dimension(:,:,:), allocatable :: constitutive_nonlocal_Cslip_66 ! elasticity matrix in Mandel notation for each instance
real(pReal), dimension(:,:,:,:,:), allocatable :: constitutive_nonlocal_Cslip_3333 ! elasticity matrix for each instance
real(pReal), dimension(:,:), allocatable :: constitutive_nonlocal_rhoEdgePos0, & ! initial edge_pos dislocation density per slip system for each family and instance
@ -145,6 +144,7 @@ integer(pInt) section, &
j, &
k, &
l, &
ns, & ! short notation for total number of active slip systems for the current instance
o, & ! index of my output
s, & ! index of my slip system
s1, & ! index of my slip system
@ -202,6 +202,7 @@ allocate(constitutive_nonlocal_nu(maxNinstance))
allocate(constitutive_nonlocal_atomicVolume(maxNinstance))
allocate(constitutive_nonlocal_D0(maxNinstance))
allocate(constitutive_nonlocal_Qsd(maxNinstance))
allocate(constitutive_nonlocal_relevantRho(maxNinstance))
allocate(constitutive_nonlocal_Cslip_66(6,6,maxNinstance))
allocate(constitutive_nonlocal_Cslip_3333(3,3,3,3,maxNinstance))
constitutive_nonlocal_CoverA = 0.0_pReal
@ -214,6 +215,7 @@ constitutive_nonlocal_Gmod = 0.0_pReal
constitutive_nonlocal_atomicVolume = 0.0_pReal
constitutive_nonlocal_D0 = 0.0_pReal
constitutive_nonlocal_Qsd = 0.0_pReal
constitutive_nonlocal_relevantRho = 0.0_pReal
constitutive_nonlocal_nu = 0.0_pReal
constitutive_nonlocal_Cslip_66 = 0.0_pReal
constitutive_nonlocal_Cslip_3333 = 0.0_pReal
@ -316,6 +318,8 @@ do
constitutive_nonlocal_D0(i) = IO_floatValue(line,positions,2)
case('qsd')
constitutive_nonlocal_Qsd(i) = IO_floatValue(line,positions,2)
case('relevantrho')
constitutive_nonlocal_relevantRho(i) = IO_floatValue(line,positions,2)
case ('interaction_slipslip')
forall (it = 1:lattice_maxNinteraction) constitutive_nonlocal_interactionSlipSlip(it,i) = IO_floatValue(line,positions,1+it)
end select
@ -327,15 +331,14 @@ enddo
constitutive_nonlocal_structure(i) = &
lattice_initializeStructure(constitutive_nonlocal_structureName(i), constitutive_nonlocal_CoverA(i)) ! our lattice structure is defined in the material.config file by the structureName (and the c/a ratio)
myStructure = constitutive_nonlocal_structure(i)
!*** sanity checks
!*** sanity checks
!*** !!! not yet complete !!!
if (constitutive_nonlocal_structure(i) < 1 .or. constitutive_nonlocal_structure(i) > 3) call IO_error(205)
if (sum(constitutive_nonlocal_Nslip(:,i)) <= 0) call IO_error(225)
if (myStructure < 1 .or. myStructure > 3) call IO_error(205)
if (sum(constitutive_nonlocal_Nslip(:,i)) <= 0_pInt) call IO_error(225)
do f = 1,lattice_maxNslipFamily
if (constitutive_nonlocal_Nslip(f,i) > 0) then
if (constitutive_nonlocal_Nslip(f,i) > 0_pInt) then
if (constitutive_nonlocal_rhoEdgePos0(f,i) < 0.0_pReal) call IO_error(220)
if (constitutive_nonlocal_rhoEdgeNeg0(f,i) < 0.0_pReal) call IO_error(220)
if (constitutive_nonlocal_rhoScrewPos0(f,i) < 0.0_pReal) call IO_error(220)
@ -343,19 +346,23 @@ enddo
if (constitutive_nonlocal_rhoEdgeDip0(f,i) < 0.0_pReal) call IO_error(220)
if (constitutive_nonlocal_rhoScrewDip0(f,i) < 0.0_pReal) call IO_error(220)
if (constitutive_nonlocal_burgersPerSlipFamily(f,i) <= 0.0_pReal) call IO_error(221)
if (constitutive_nonlocal_v0PerSlipFamily(f,i) <= 0.0_pReal) call IO_error(-1)
if (constitutive_nonlocal_lambda0PerSlipFamily(f,i) <= 0.0_pReal) call IO_error(-1)
if (constitutive_nonlocal_dDipMinEdgePerSlipFamily(f,i) <= 0.0_pReal) call IO_error(-1)
if (constitutive_nonlocal_dDipMinScrewPerSlipFamily(f,i) <= 0.0_pReal) call IO_error(-1)
if (constitutive_nonlocal_v0PerSlipFamily(f,i) <= 0.0_pReal) call IO_error(226)
if (constitutive_nonlocal_lambda0PerSlipFamily(f,i) <= 0.0_pReal) call IO_error(227)
if (constitutive_nonlocal_dDipMinEdgePerSlipFamily(f,i) <= 0.0_pReal) call IO_error(228)
if (constitutive_nonlocal_dDipMinScrewPerSlipFamily(f,i) <= 0.0_pReal) call IO_error(228)
endif
enddo
if (sum(constitutive_nonlocal_interactionSlipSlip(:,i)) <= 0.0_pReal) call IO_error(-1)
if (any(constitutive_nonlocal_interactionSlipSlip(1:maxval(lattice_interactionSlipSlip(:,:,myStructure)),i) < 1.0_pReal)) &
call IO_error(229)
if (constitutive_nonlocal_atomicVolume(i) <= 0.0_pReal) call IO_error(230)
if (constitutive_nonlocal_D0(i) <= 0.0_pReal) call IO_error(231)
if (constitutive_nonlocal_Qsd(i) <= 0.0_pReal) call IO_error(232)
if (constitutive_nonlocal_relevantRho(i) <= 0.0_pReal) call IO_error(233)
!*** determine total number of active slip systems
constitutive_nonlocal_Nslip(:,i) &
= min( lattice_NslipSystem(:, constitutive_nonlocal_structure(i)), constitutive_nonlocal_Nslip(:,i) ) ! we can't use more slip systems per family than specified in lattice
constitutive_nonlocal_Nslip(:,i) = min( lattice_NslipSystem(:, myStructure), constitutive_nonlocal_Nslip(:,i) ) ! we can't use more slip systems per family than specified in lattice
constitutive_nonlocal_totalNslip(i) = sum(constitutive_nonlocal_Nslip(:,i))
enddo
@ -407,8 +414,10 @@ do i = 1,maxNinstance
!*** determine size of state array
constitutive_nonlocal_sizeState(i) = (size(constitutive_nonlocal_stateList)-1) * constitutive_nonlocal_totalNslip(i) + 6_pInt ! the size of the list of states times the number of active slip systems gives the required size for the state array
constitutive_nonlocal_sizeDotState(i) = size(constitutive_nonlocal_stateListBasic) * constitutive_nonlocal_totalNslip(i) ! the size of the list of basic states times the number of active slip systems gives the required size for the dotState array
ns = constitutive_nonlocal_totalNslip(i)
constitutive_nonlocal_sizeState(i) = size(constitutive_nonlocal_listBasicStates) * ns &
+ ( size(constitutive_nonlocal_listDependentStates) - 1_pInt ) * ns + 6_pInt
constitutive_nonlocal_sizeDotState(i) = size(constitutive_nonlocal_listBasicStates) * ns
!*** determine size of postResults array
@ -528,7 +537,6 @@ pure function constitutive_nonlocal_stateInit(myInstance)
use prec, only: pReal, &
pInt
use lattice, only: lattice_maxNslipFamily
use IO, only: IO_error
implicit none
!*** input variables
@ -610,6 +618,30 @@ endfunction
!*********************************************************************
!* relevant microstructural state *
!*********************************************************************
pure function constitutive_nonlocal_relevantState(myInstance)
use prec, only: pReal, &
pInt
implicit none
!*** input variables
integer(pInt), intent(in) :: myInstance ! number specifying the current instance of the constitution
!*** output variables
real(pReal), dimension(constitutive_nonlocal_sizeState(myInstance)) :: &
constitutive_nonlocal_relevantState ! relevant state values for the current instance of this constitution
!*** local variables
constitutive_nonlocal_relevantState = constitutive_nonlocal_relevantRho(myInstance)
endfunction
!*********************************************************************
!* calculates homogenized elacticity matrix *
!*********************************************************************
@ -838,7 +870,7 @@ do n = 1,FE_NipNeighbors(mesh_element(2,el))
sigma(3,1) = 0.0_pReal
! coordinate transformation from the slip coordinate system to the lattice coordinate system
Tdislocation_v = Tdislocation_v + math_Mandel33to6( math_mul33x33(transpose(transform), math_mul33x33(sigma, transform) ) )
Tdislocation_v = Tdislocation_v + math_Mandel33to6( math_mul33x33(transform, math_mul33x33(sigma, transpose(transform)) ) )
! if (debugger) write(6,'(a15,3(i3,x),/,3(3(f12.3,x)/))') 'sigma / MPa at ',g,ip,el, sigma/1e6
! if (debugger) write(6,'(a15,3(i3,x),/,3(3(f12.3,x)/))') 'Tdislocation / MPa at ',g,ip,el, math_Mandel6to33(Tdislocation_v/1e6)
enddo
@ -1142,7 +1174,7 @@ if (debugger) then
write(6,*)
write(6,'(a,/,12(f12.5,x),/)') 'tauSlip / MPa', tauSlip/1e6_pReal
write(6,'(a,/,12(f12.5,x),/)') 'tauSlipThreshold / MPa', tauSlipThreshold/1e6_pReal
! write(6,'(a,/,12(e10.3,x),/)') 'v', v
write(6,'(a,/,12(e10.3,x),/)') 'v', v
write(6,'(a,/,4(12(f12.5,x),/))') 'gdot / 1e-3', gdot*1e3_pReal
write(6,'(a,/,(12(f12.5,x),/))') 'gdot total/ 1e-3', sum(gdot,2)*1e3_pReal
write(6,'(a,/,6(12(e12.5,x),/))') 'dislocation multiplication', &
@ -1221,7 +1253,7 @@ dDipMaxDot(:,2) = (dDipMax(:,2) - dDipMax0(:,2)) / subdt
forall (c=1:2) &
rhoDotTransfer(:,c) = 2.0_pReal * dDipMax(:,c) / constitutive_nonlocal_burgersPerSlipSystem(:,myInstance) &
* ( rho(:,2*c-1)*gdot(:,2*c) + rho(:,2*c)*gdot(:,2*c-1) )
* ( rho(:,2*c-1)*abs(gdot(:,2*c)) + rho(:,2*c)*abs(gdot(:,2*c-1)) )
if (debugger) write(6,'(a,/,6(12(e12.5,x),/))') 'dipole formation by glide', &
-rhoDotTransfer*subdt,-rhoDotTransfer*subdt,2.0_pReal*rhoDotTransfer*subdt
@ -1234,7 +1266,7 @@ rhoDipDot = rhoDipDot + 2.0_pReal * rhoDotTransfer
forall (c=1:2) &
rhoDotTransfer(:,c) = - 2.0_pReal * dDipMin(:,c) / constitutive_nonlocal_burgersPerSlipSystem(:,myInstance) &
* ( rho(:,2*c-1)*gdot(:,2*c) + rho(:,2*c)*gdot(:,2*c-1) )
* ( rho(:,2*c-1)*abs(gdot(:,2*c)) + rho(:,2*c)*abs(gdot(:,2*c-1)) )
if (debugger) write(6,'(a,/,6(12(e12.5,x),/))') 'athermal dipole annihilation', &
0.0_pReal*rhoDotTransfer,0.0_pReal*rhoDotTransfer,2.0_pReal*rhoDotTransfer*subdt

View File

@ -46,6 +46,7 @@
! interaction_sliptwin 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16
! interaction_twinslip 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16
! interaction_twintwin 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20
! relevantResistance 1e2
MODULE constitutive_phenopowerlaw
@ -107,6 +108,8 @@ MODULE constitutive_phenopowerlaw
real(pReal), dimension(:), allocatable :: constitutive_phenopowerlaw_w0_slip
real(pReal), dimension(:), allocatable :: constitutive_phenopowerlaw_relevantResistance
CONTAINS
!****************************************
!* - constitutive_init
@ -205,7 +208,11 @@ subroutine constitutive_phenopowerlaw_init(file)
constitutive_phenopowerlaw_interaction_twinslip = 0.0_pReal
constitutive_phenopowerlaw_interaction_twintwin = 0.0_pReal
allocate(constitutive_phenopowerlaw_w0_slip(maxNinstance)) ; constitutive_phenopowerlaw_w0_slip = 0.0_pReal
allocate(constitutive_phenopowerlaw_w0_slip(maxNinstance))
constitutive_phenopowerlaw_w0_slip = 0.0_pReal
allocate(constitutive_phenopowerlaw_relevantResistance(maxNinstance))
constitutive_phenopowerlaw_relevantResistance = 0.0_pReal
rewind(file)
line = ''
@ -282,6 +289,8 @@ subroutine constitutive_phenopowerlaw_init(file)
constitutive_phenopowerlaw_h0_twinslip(i) = IO_floatValue(line,positions,2)
case ('h0_twintwin')
constitutive_phenopowerlaw_h0_twintwin(i) = IO_floatValue(line,positions,2)
case ('relevantresistance')
constitutive_phenopowerlaw_relevantResistance(i) = IO_floatValue(line,positions,2)
case ('interaction_slipslip')
forall (j = 1:lattice_maxNinteraction) &
constitutive_phenopowerlaw_interaction_slipslip(j,i) = IO_floatValue(line,positions,1+j)
@ -308,21 +317,22 @@ subroutine constitutive_phenopowerlaw_init(file)
constitutive_phenopowerlaw_totalNslip(i) = sum(constitutive_phenopowerlaw_Nslip(:,i)) ! how many slip systems altogether
constitutive_phenopowerlaw_totalNtwin(i) = sum(constitutive_phenopowerlaw_Ntwin(:,i)) ! how many twin systems altogether
if (constitutive_phenopowerlaw_structure(i) < 1 ) call IO_error(205)
if (constitutive_phenopowerlaw_structure(i) < 1 ) call IO_error(205)
if (any(constitutive_phenopowerlaw_tau0_slip(:,i) < 0.0_pReal .and. &
constitutive_phenopowerlaw_Nslip(:,i) > 0)) call IO_error(210)
if (constitutive_phenopowerlaw_gdot0_slip(i) <= 0.0_pReal) call IO_error(211)
if (constitutive_phenopowerlaw_n_slip(i) <= 0.0_pReal) call IO_error(212)
constitutive_phenopowerlaw_Nslip(:,i) > 0)) call IO_error(210)
if (constitutive_phenopowerlaw_gdot0_slip(i) <= 0.0_pReal) call IO_error(211)
if (constitutive_phenopowerlaw_n_slip(i) <= 0.0_pReal) call IO_error(212)
if (any(constitutive_phenopowerlaw_tausat_slip(:,i) <= 0.0_pReal .and. &
constitutive_phenopowerlaw_Nslip(:,i) > 0)) call IO_error(213)
constitutive_phenopowerlaw_Nslip(:,i) > 0)) call IO_error(213)
if (any(constitutive_phenopowerlaw_w0_slip(i) == 0.0_pReal .and. &
constitutive_phenopowerlaw_Nslip(:,i) > 0)) call IO_error(214)
constitutive_phenopowerlaw_Nslip(:,i) > 0)) call IO_error(214)
if (any(constitutive_phenopowerlaw_tau0_twin(:,i) < 0.0_pReal .and. &
constitutive_phenopowerlaw_Ntwin(:,i) > 0)) call IO_error(210)
if ( constitutive_phenopowerlaw_gdot0_twin(i) <= 0.0_pReal .and. &
any(constitutive_phenopowerlaw_Ntwin(:,i) > 0)) call IO_error(211)
constitutive_phenopowerlaw_Ntwin(:,i) > 0)) call IO_error(210)
if ( constitutive_phenopowerlaw_gdot0_twin(i) <= 0.0_pReal .and. &
any(constitutive_phenopowerlaw_Ntwin(:,i) > 0)) call IO_error(211)
if ( constitutive_phenopowerlaw_n_twin(i) <= 0.0_pReal .and. &
any(constitutive_phenopowerlaw_Ntwin(:,i) > 0)) call IO_error(212)
any(constitutive_phenopowerlaw_Ntwin(:,i) > 0)) call IO_error(212)
if (constitutive_phenopowerlaw_relevantResistance(i) <= 0.0_pReal) call IO_error(242)
enddo
@ -492,6 +502,29 @@ function constitutive_phenopowerlaw_stateInit(myInstance)
endfunction
!*********************************************************************
!* relevant microstructural state *
!*********************************************************************
pure function constitutive_phenopowerlaw_relevantState(myInstance)
use prec, only: pReal, &
pInt
implicit none
!*** input variables
integer(pInt), intent(in) :: myInstance ! number specifying the current instance of the constitution
!*** output variables
real(pReal), dimension(constitutive_phenopowerlaw_sizeState(myInstance)) :: &
constitutive_phenopowerlaw_relevantState ! relevant state values for the current instance of this constitution
!*** local variables
constitutive_phenopowerlaw_relevantState = constitutive_phenopowerlaw_relevantResistance(myInstance)
endfunction
function constitutive_phenopowerlaw_homogenizedC(state,ipc,ip,el)
!*********************************************************************
!* homogenized elacticity matrix *

View File

@ -701,6 +701,7 @@ endsubroutine
constitutive_sizeDotState, &
constitutive_subState0, &
constitutive_state, &
constitutive_relevantState, &
constitutive_microstructure
use debug, only: debugger, &
debug_cumDotStateCalls, &
@ -751,7 +752,7 @@ endsubroutine
! setting flag to true if state is below relative tolerance, otherwise set it to false
crystallite_updateState = all( constitutive_state(g,i,e)%p(1:mySize) == 0.0_pReal &
crystallite_updateState = all( constitutive_state(g,i,e)%p(1:mySize) < constitutive_relevantState(g,i,e)%p(1:mySize) &
.or. abs(residuum) < rTol_crystalliteState*abs(constitutive_state(g,i,e)%p(1:mySize)))
if (debugger) then
!$OMP CRITICAL (write2out)

View File

@ -58,7 +58,8 @@ gdot0 0.001
n 20
h0 75e6
tausat 63e6
w0 1
w0 2.25
relevantResistance 1
[Aluminum_phenopowerlaw]
# slip only
@ -102,6 +103,7 @@ interaction_slipslip 1 1 1.4 1.4 1.4 1.4
interaction_sliptwin 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1
interaction_twinslip 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1
interaction_twintwin 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1
relevantResistance 1
[Aluminum_nonlocal]
@ -119,26 +121,27 @@ constitution nonlocal
(output) resistance
lattice_structure fcc
Nslip 12 0 0 0 # per family
Nslip 12 0 0 0 # number of slip systems per family
c11 106.75e9
c12 60.41e9
c44 28.34e9
burgers 2.86e-10 0 0 0 # Burgers vector in m
rhoEdgePos0 1.0e10 0 0 0 # Initial positive edge dislocation density in m/m**3
rhoEdgeNeg0 1.0e10 0 0 0 # Initial negative edge dislocation density in m/m**3
rhoScrewPos0 1.0e10 0 0 0 # Initial positive screw dislocation density in m/m**3
rhoScrewNeg0 1.0e10 0 0 0 # Initial negative screw dislocation density in m/m**3
rhoEdgeDip0 0 0 0 0 # Initial edge dipole dislocation density in m/m**3
rhoScrewDip0 0 0 0 0 # Initial screw dipole dislocation density in m/m**3
rhoEdgePos0 1e10 0 0 0 # Initial positive edge dislocation density in m/m**3
rhoEdgeNeg0 1e10 0 0 0 # Initial negative edge dislocation density in m/m**3
rhoScrewPos0 1e10 0 0 0 # Initial positive screw dislocation density in m/m**3
rhoScrewNeg0 1e10 0 0 0 # Initial negative screw dislocation density in m/m**3
rhoEdgeDip0 1 0 0 0 # Initial edge dipole dislocation density in m/m**3
rhoScrewDip0 1 0 0 0 # Initial screw dipole dislocation density in m/m**3
v0 1e-4 0 0 0 # prefactor for dislocation velocity
dDipMinEdge 1e-9 0 0 0 # minimum distance for stable edge dipoles in m
dDipMinScrew 1e-9 0 0 0 # minimum distance for stable screw dipoles in m
lambda0 100 0 0 0 # prefactor for mean free path
atomicVolume 1.7e-29
D0 1e-4
Qsd 2.3e-19
lambda0 100 0 0 0 # prefactor for mean free path
D0 1e-4 # prefactor for self-diffusion coefficient
Qsd 2.3e-19 # activation energy for dislocation climb
relevantRho 1e3 # dislocation density considered relevant
interaction_SlipSlip 1.0 2.2 3.0 1.6 3.8 4.5 # Dislocation interaction coefficient
[BCC_Ferrite]
@ -170,6 +173,7 @@ interaction_sliptwin 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1
interaction_twinslip 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1
interaction_twintwin 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1
w0_slip 1.0
relevantResistance 1
[BCC_Martensite]
constitution phenopowerlaw
@ -200,6 +204,7 @@ interaction_sliptwin 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1
interaction_twinslip 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1
interaction_twintwin 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1
w0_slip 4.0
relevantResistance 1
#####################
<texture>