Modification of basic states for phase transformation
This commit is contained in:
parent
591a82de26
commit
1da874a6db
|
@ -30,7 +30,7 @@ module constitutive_dislotwin
|
|||
CONSTITUTIVE_DISLOTWIN_listBasicTwinStates = &
|
||||
['twinFraction', 'accsheartwin']
|
||||
|
||||
character(len=12), dimension(1), parameter, private :: &
|
||||
character(len=13), dimension(1), parameter, private :: &
|
||||
CONSTITUTIVE_DISLOTWIN_listBasicTransStates = &
|
||||
['transFraction']
|
||||
|
||||
|
@ -953,12 +953,12 @@ subroutine constitutive_dislotwin_stateInit(ph,instance)
|
|||
forall (i = 1_pInt:ns) &
|
||||
invLambdaSlip0(i) = sqrt(dot_product((rhoEdge0+rhoEdgeDip0),constitutive_dislotwin_forestProjectionEdge(1:ns,i,instance)))/ &
|
||||
constitutive_dislotwin_CLambdaSlipPerSlipSystem(i,instance)
|
||||
tempState(3_pInt*ns+2_pInt*nt+1:4_pInt*ns+2_pInt*nt) = invLambdaSlip0
|
||||
tempState(3_pInt*ns+2_pInt*nt+nr+1:4_pInt*ns+2_pInt*nt+nr) = invLambdaSlip0
|
||||
|
||||
forall (i = 1_pInt:ns) &
|
||||
MeanFreePathSlip0(i) = &
|
||||
constitutive_dislotwin_GrainSize(instance)/(1.0_pReal+invLambdaSlip0(i)*constitutive_dislotwin_GrainSize(instance))
|
||||
tempState(5_pInt*ns+3_pInt*nt+1:6_pInt*ns+3_pInt*nt) = MeanFreePathSlip0
|
||||
tempState(5_pInt*ns+3_pInt*nt+nr+1:6_pInt*ns+3_pInt*nt+nr) = MeanFreePathSlip0
|
||||
|
||||
forall (i = 1_pInt:ns) &
|
||||
tauSlipThreshold0(i) = &
|
||||
|
@ -967,8 +967,6 @@ subroutine constitutive_dislotwin_stateInit(ph,instance)
|
|||
|
||||
tempState(6_pInt*ns+4_pInt*nt+nr+1:7_pInt*ns+4_pInt*nt+nr) = tauSlipThreshold0
|
||||
|
||||
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
! initialize dependent twin microstructural variables
|
||||
forall (j = 1_pInt:nt) &
|
||||
|
@ -1123,13 +1121,14 @@ subroutine constitutive_dislotwin_microstructure(temperature,ipc,ip,el)
|
|||
ns = constitutive_dislotwin_totalNslip(instance)
|
||||
nt = constitutive_dislotwin_totalNtwin(instance)
|
||||
nr = constitutive_dislotwin_totalNtrans(instance)
|
||||
|
||||
!BASIC STATES
|
||||
!* State: 1 : ns rho_edge
|
||||
!* State: ns+1 : 2*ns rho_dipole
|
||||
!* State: 2*ns+1 : 3*ns accumulated shear due to slip
|
||||
!* State: 3*ns+1 : 3*ns+nt f
|
||||
!* State: 3*ns+nt+1 : 3*ns+2*nt accumulated shear due to twin
|
||||
!* State: 3*ns+2*nt+1 : 3*ns+2*nt+nr transformed volume fraction
|
||||
!* State: 1 : ns rho_edge
|
||||
!* State: ns+1 : 2*ns rho_dipole
|
||||
!* State: 2*ns+1 : 3*ns accumulated shear due to slip
|
||||
!* State: 3*ns+1 : 3*ns+nt f
|
||||
!* State: 3*ns+nt+1 : 3*ns+2*nt accumulated shear due to twin
|
||||
!* State: 3*ns+2*nt+1 : 3*ns+2*nt+nr transformed volume fraction
|
||||
!DEPENDENT STATES
|
||||
!* State: 3*ns+2*nt+nr+1 : 4*ns+2*nt+nr 1/lambda_slip
|
||||
!* State: 4*ns+2*nt+nr+1 : 5*ns+2*nt+nr 1/lambda_sliptwin
|
||||
|
@ -1158,6 +1157,7 @@ subroutine constitutive_dislotwin_microstructure(temperature,ipc,ip,el)
|
|||
sqrt(dot_product((plasticState(ph)%state(1:ns,of)+plasticState(ph)%state(ns+1_pInt:2_pInt*ns,of)),&
|
||||
constitutive_dislotwin_forestProjectionEdge(1:ns,s,instance)))/ &
|
||||
constitutive_dislotwin_CLambdaSlipPerSlipSystem(s,instance)
|
||||
|
||||
!* 1/mean free distance between 2 twin stacks from different systems seen by a moving dislocation
|
||||
!$OMP CRITICAL (evilmatmul)
|
||||
plasticState(ph)%state((4_pInt*ns+2_pInt*nt+nr+1_pInt):(5_pInt*ns+2_pInt*nt+nr), of) = 0.0_pReal
|
||||
|
@ -1255,8 +1255,10 @@ subroutine constitutive_dislotwin_LpAndItsTangent(Lp,dLp_dTstar,Tstar_v,Temperat
|
|||
lattice_Stwin_v, &
|
||||
lattice_maxNslipFamily,&
|
||||
lattice_maxNtwinFamily, &
|
||||
lattice_maxNtransFamily, &
|
||||
lattice_NslipSystem, &
|
||||
lattice_NtwinSystem, &
|
||||
lattice_NtransSystem, &
|
||||
lattice_shearTwin, &
|
||||
lattice_structure, &
|
||||
lattice_fcc_twinNucleationSlipPair, &
|
||||
|
@ -1269,7 +1271,7 @@ subroutine constitutive_dislotwin_LpAndItsTangent(Lp,dLp_dTstar,Tstar_v,Temperat
|
|||
real(pReal), dimension(3,3), intent(out) :: Lp
|
||||
real(pReal), dimension(9,9), intent(out) :: dLp_dTstar
|
||||
|
||||
integer(pInt) :: instance,ph,of,ns,nt,f,i,j,k,l,m,n,index_myFamily,s1,s2
|
||||
integer(pInt) :: instance,ph,of,ns,nt,nr,f,i,j,k,l,m,n,index_myFamily,s1,s2
|
||||
real(pReal) :: sumf,StressRatio_p,StressRatio_pminus1,StressRatio_r,BoltzmannRatio,DotGamma0,Ndot0
|
||||
real(pReal), dimension(3,3,3,3) :: dLp_dTstar3333
|
||||
real(pReal), dimension(constitutive_dislotwin_totalNslip(phase_plasticityInstance(material_phase(ipc,ip,el)))) :: &
|
||||
|
@ -1305,6 +1307,7 @@ subroutine constitutive_dislotwin_LpAndItsTangent(Lp,dLp_dTstar,Tstar_v,Temperat
|
|||
instance = phase_plasticityInstance(ph)
|
||||
ns = constitutive_dislotwin_totalNslip(instance)
|
||||
nt = constitutive_dislotwin_totalNtwin(instance)
|
||||
nr = constitutive_dislotwin_totalNtrans(instance)
|
||||
|
||||
!* Total twin volume fraction
|
||||
sumf = sum(plasticState(ph)%state((3_pInt*ns+1_pInt):(3_pInt*ns+nt), of)) ! safe for nt == 0
|
||||
|
@ -1326,12 +1329,12 @@ subroutine constitutive_dislotwin_LpAndItsTangent(Lp,dLp_dTstar,Tstar_v,Temperat
|
|||
!* Resolved shear stress on slip system
|
||||
tau_slip(j) = dot_product(Tstar_v,lattice_Sslip_v(:,1,index_myFamily+i,ph))
|
||||
|
||||
if((abs(tau_slip(j))-plasticState(ph)%state(6*ns+4*nt+j, of)) > tol_math_check) then
|
||||
if((abs(tau_slip(j))-plasticState(ph)%state(6*ns+4*nt+nr+j, of)) > tol_math_check) then
|
||||
!* Stress ratios
|
||||
StressRatio_p = ((abs(tau_slip(j))- plasticState(ph)%state(6*ns+4*nt+j, of))/&
|
||||
StressRatio_p = ((abs(tau_slip(j))- plasticState(ph)%state(6*ns+4*nt+nr+j, of))/&
|
||||
(constitutive_dislotwin_SolidSolutionStrength(instance)+constitutive_dislotwin_tau_peierlsPerSlipFamily(f,instance)))&
|
||||
**constitutive_dislotwin_pPerSlipFamily(f,instance)
|
||||
StressRatio_pminus1 = ((abs(tau_slip(j))-plasticState(ph)%state(6*ns+4*nt+j, of))/&
|
||||
StressRatio_pminus1 = ((abs(tau_slip(j))-plasticState(ph)%state(6*ns+4*nt+nr+j, of))/&
|
||||
(constitutive_dislotwin_SolidSolutionStrength(instance)+constitutive_dislotwin_tau_peierlsPerSlipFamily(f,instance)))&
|
||||
**(constitutive_dislotwin_pPerSlipFamily(f,instance)-1.0_pReal)
|
||||
!* Boltzmann ratio
|
||||
|
@ -1437,7 +1440,7 @@ subroutine constitutive_dislotwin_LpAndItsTangent(Lp,dLp_dTstar,Tstar_v,Temperat
|
|||
|
||||
!* Stress ratios
|
||||
if (tau_twin(j) > tol_math_check) then
|
||||
StressRatio_r = (plasticState(ph)%state(7*ns+4*nt+j, of)/tau_twin(j))**constitutive_dislotwin_rPerTwinFamily(f,instance)
|
||||
StressRatio_r = (plasticState(ph)%state(7*ns+4*nt+nr+j, of)/tau_twin(j))**constitutive_dislotwin_rPerTwinFamily(f,instance)
|
||||
!* Shear rates and their derivatives due to twin
|
||||
select case(lattice_structure(ph))
|
||||
case (LATTICE_fcc_ID)
|
||||
|
@ -1457,7 +1460,7 @@ subroutine constitutive_dislotwin_LpAndItsTangent(Lp,dLp_dTstar,Tstar_v,Temperat
|
|||
end select
|
||||
gdot_twin(j) = &
|
||||
(constitutive_dislotwin_MaxTwinFraction(instance)-sumf)*lattice_shearTwin(index_myFamily+i,ph)*&
|
||||
plasticState(ph)%state(7*ns+5*nt+j, of)*Ndot0*exp(-StressRatio_r)
|
||||
plasticState(ph)%state(7*ns+5*nt+nr+j, of)*Ndot0*exp(-StressRatio_r)
|
||||
dgdot_dtautwin(j) = ((gdot_twin(j)*constitutive_dislotwin_rPerTwinFamily(f,instance))/tau_twin(j))*StressRatio_r
|
||||
endif
|
||||
|
||||
|
@ -1501,8 +1504,10 @@ subroutine constitutive_dislotwin_dotState(Tstar_v,Temperature,ipc,ip,el)
|
|||
lattice_Stwin_v, &
|
||||
lattice_maxNslipFamily, &
|
||||
lattice_maxNtwinFamily, &
|
||||
lattice_maxNtransFamily, &
|
||||
lattice_NslipSystem, &
|
||||
lattice_NtwinSystem, &
|
||||
lattice_NtransSystem, &
|
||||
lattice_sheartwin, &
|
||||
lattice_mu, &
|
||||
lattice_structure, &
|
||||
|
@ -1520,7 +1525,7 @@ subroutine constitutive_dislotwin_dotState(Tstar_v,Temperature,ipc,ip,el)
|
|||
ip, & !< integration point
|
||||
el !< element
|
||||
|
||||
integer(pInt) :: instance,ns,nt,f,i,j,index_myFamily,s1,s2, &
|
||||
integer(pInt) :: instance,ns,nt,nr,f,i,j,index_myFamily,s1,s2, &
|
||||
ph, &
|
||||
of
|
||||
real(pReal) :: sumf,StressRatio_p,StressRatio_pminus1,BoltzmannRatio,DotGamma0,&
|
||||
|
@ -1538,7 +1543,7 @@ subroutine constitutive_dislotwin_dotState(Tstar_v,Temperature,ipc,ip,el)
|
|||
instance = phase_plasticityInstance(ph)
|
||||
ns = constitutive_dislotwin_totalNslip(instance)
|
||||
nt = constitutive_dislotwin_totalNtwin(instance)
|
||||
|
||||
nr = constitutive_dislotwin_totalNtrans(instance)
|
||||
|
||||
!* Total twin volume fraction
|
||||
sumf = sum(plasticState(ph)%state((3_pInt*ns+1_pInt):(3_pInt*ns+nt), of)) ! safe for nt == 0
|
||||
|
@ -1556,12 +1561,12 @@ subroutine constitutive_dislotwin_dotState(Tstar_v,Temperature,ipc,ip,el)
|
|||
!* Resolved shear stress on slip system
|
||||
tau_slip(j) = dot_product(Tstar_v,lattice_Sslip_v(:,1,index_myFamily+i,ph))
|
||||
|
||||
if((abs(tau_slip(j))-plasticState(ph)%state(6*ns+4*nt+j, of)) > tol_math_check) then
|
||||
if((abs(tau_slip(j))-plasticState(ph)%state(6*ns+4*nt+nr+j, of)) > tol_math_check) then
|
||||
!* Stress ratios
|
||||
StressRatio_p = ((abs(tau_slip(j))-plasticState(ph)%state(6*ns+4*nt+j, of))/&
|
||||
StressRatio_p = ((abs(tau_slip(j))-plasticState(ph)%state(6*ns+4*nt+nr+j, of))/&
|
||||
(constitutive_dislotwin_SolidSolutionStrength(instance)+constitutive_dislotwin_tau_peierlsPerSlipFamily(f,instance)))&
|
||||
**constitutive_dislotwin_pPerSlipFamily(f,instance)
|
||||
StressRatio_pminus1 = ((abs(tau_slip(j))-plasticState(ph)%state(6*ns+4*nt+j, of))/&
|
||||
StressRatio_pminus1 = ((abs(tau_slip(j))-plasticState(ph)%state(6*ns+4*nt+nr+j, of))/&
|
||||
(constitutive_dislotwin_SolidSolutionStrength(instance)+constitutive_dislotwin_tau_peierlsPerSlipFamily(f,instance)))&
|
||||
**(constitutive_dislotwin_pPerSlipFamily(f,instance)-1.0_pReal)
|
||||
!* Boltzmann ratio
|
||||
|
@ -1578,7 +1583,7 @@ subroutine constitutive_dislotwin_dotState(Tstar_v,Temperature,ipc,ip,el)
|
|||
!* Multiplication
|
||||
DotRhoMultiplication(j) = abs(gdot_slip(j))/&
|
||||
(constitutive_dislotwin_burgersPerSlipSystem(j,instance)* &
|
||||
plasticState(ph)%state(5*ns+3*nt+j, of))
|
||||
plasticState(ph)%state(5*ns+3*nt+nr+j, of))
|
||||
!* Dipole formation
|
||||
EdgeDipMinDistance = &
|
||||
constitutive_dislotwin_CEdgeDipMinDistance(instance)*constitutive_dislotwin_burgersPerSlipSystem(j,instance)
|
||||
|
@ -1588,7 +1593,8 @@ subroutine constitutive_dislotwin_dotState(Tstar_v,Temperature,ipc,ip,el)
|
|||
EdgeDipDistance(j) = &
|
||||
(3.0_pReal*lattice_mu(ph)*constitutive_dislotwin_burgersPerSlipSystem(j,instance))/&
|
||||
(16.0_pReal*pi*abs(tau_slip(j)))
|
||||
if (EdgeDipDistance(j)>plasticState(ph)%state(5*ns+3*nt+j, of)) EdgeDipDistance(j)=plasticState(ph)%state(5*ns+3*nt+j, of)
|
||||
if (EdgeDipDistance(j)>plasticState(ph)%state(5*ns+3*nt+nr+j, of)) &
|
||||
EdgeDipDistance(j)=plasticState(ph)%state(5*ns+3*nt+nr+j, of)
|
||||
if (EdgeDipDistance(j)<EdgeDipMinDistance) EdgeDipDistance(j)=EdgeDipMinDistance
|
||||
DotRhoDipFormation(j) = &
|
||||
((2.0_pReal*EdgeDipDistance(j))/constitutive_dislotwin_burgersPerSlipSystem(j,instance))*&
|
||||
|
@ -1645,7 +1651,7 @@ subroutine constitutive_dislotwin_dotState(Tstar_v,Temperature,ipc,ip,el)
|
|||
tau_twin(j) = dot_product(Tstar_v,lattice_Stwin_v(:,index_myFamily+i,ph))
|
||||
!* Stress ratios
|
||||
if (tau_twin(j) > tol_math_check) then
|
||||
StressRatio_r = (plasticState(ph)%state(7*ns+4*nt+j, of)/tau_twin(j))**constitutive_dislotwin_rPerTwinFamily(f,instance)
|
||||
StressRatio_r = (plasticState(ph)%state(7*ns+4*nt+nr+j, of)/tau_twin(j))**constitutive_dislotwin_rPerTwinFamily(f,instance)
|
||||
!* Shear rates and their derivatives due to twin
|
||||
|
||||
select case(lattice_structure(ph))
|
||||
|
@ -1666,7 +1672,7 @@ subroutine constitutive_dislotwin_dotState(Tstar_v,Temperature,ipc,ip,el)
|
|||
end select
|
||||
plasticState(ph)%dotState(3_pInt*ns+j, of) = &
|
||||
(constitutive_dislotwin_MaxTwinFraction(instance)-sumf)*&
|
||||
plasticState(ph)%state(7_pInt*ns+5_pInt*nt+j, of)*Ndot0*exp(-StressRatio_r)
|
||||
plasticState(ph)%state(7_pInt*ns+5_pInt*nt+nr+j, of)*Ndot0*exp(-StressRatio_r)
|
||||
!* Dotstate for accumulated shear due to twin
|
||||
plasticState(ph)%dotState(3_pInt*ns+nt+j, of) = plasticState(ph)%dotState(3_pInt*ns+j, of) * &
|
||||
lattice_sheartwin(index_myfamily+i,ph)
|
||||
|
@ -1703,8 +1709,10 @@ function constitutive_dislotwin_postResults(Tstar_v,Temperature,ipc,ip,el)
|
|||
lattice_Stwin_v, &
|
||||
lattice_maxNslipFamily, &
|
||||
lattice_maxNtwinFamily, &
|
||||
lattice_maxNtransFamily, &
|
||||
lattice_NslipSystem, &
|
||||
lattice_NtwinSystem, &
|
||||
lattice_NtransSystem, &
|
||||
lattice_shearTwin, &
|
||||
lattice_mu, &
|
||||
lattice_structure, &
|
||||
|
@ -1725,7 +1733,7 @@ function constitutive_dislotwin_postResults(Tstar_v,Temperature,ipc,ip,el)
|
|||
constitutive_dislotwin_postResults
|
||||
integer(pInt) :: &
|
||||
instance,phase,&
|
||||
ns,nt,&
|
||||
ns,nt,nr,&
|
||||
f,o,i,c,j,index_myFamily,&
|
||||
s1,s2, &
|
||||
ph, &
|
||||
|
@ -1744,7 +1752,7 @@ function constitutive_dislotwin_postResults(Tstar_v,Temperature,ipc,ip,el)
|
|||
instance = phase_plasticityInstance(ph)
|
||||
ns = constitutive_dislotwin_totalNslip(instance)
|
||||
nt = constitutive_dislotwin_totalNtwin(instance)
|
||||
|
||||
nr = constitutive_dislotwin_totalNtrans(instance)
|
||||
|
||||
!* Total twin volume fraction
|
||||
sumf = sum(plasticState(ph)%state((3_pInt*ns+1_pInt):(3_pInt*ns+nt), of)) ! safe for nt == 0
|
||||
|
@ -1775,13 +1783,13 @@ function constitutive_dislotwin_postResults(Tstar_v,Temperature,ipc,ip,el)
|
|||
!* Resolved shear stress on slip system
|
||||
tau = dot_product(Tstar_v,lattice_Sslip_v(:,1,index_myFamily+i,ph))
|
||||
!* Stress ratios
|
||||
if((abs(tau)-plasticState(ph)%state(6*ns+4*nt+j, of)) > tol_math_check) then
|
||||
if((abs(tau)-plasticState(ph)%state(6*ns+4*nt+nr+j, of)) > tol_math_check) then
|
||||
!* Stress ratios
|
||||
StressRatio_p = ((abs(tau)-plasticState(ph)%state(6*ns+4*nt+j, of))/&
|
||||
StressRatio_p = ((abs(tau)-plasticState(ph)%state(6*ns+4*nt+nr+j, of))/&
|
||||
(constitutive_dislotwin_SolidSolutionStrength(instance)+&
|
||||
constitutive_dislotwin_tau_peierlsPerSlipFamily(f,instance)))&
|
||||
**constitutive_dislotwin_pPerSlipFamily(f,instance)
|
||||
StressRatio_pminus1 = ((abs(tau)-plasticState(ph)%state(6*ns+4*nt+j, of))/&
|
||||
StressRatio_pminus1 = ((abs(tau)-plasticState(ph)%state(6*ns+4*nt+nr+j, of))/&
|
||||
(constitutive_dislotwin_SolidSolutionStrength(instance)+&
|
||||
constitutive_dislotwin_tau_peierlsPerSlipFamily(f,instance)))&
|
||||
**(constitutive_dislotwin_pPerSlipFamily(f,instance)-1.0_pReal)
|
||||
|
@ -1808,7 +1816,7 @@ function constitutive_dislotwin_postResults(Tstar_v,Temperature,ipc,ip,el)
|
|||
c = c + ns
|
||||
case (mfp_slip_ID)
|
||||
constitutive_dislotwin_postResults(c+1_pInt:c+ns) =&
|
||||
plasticState(ph)%state((5_pInt*ns+3_pInt*nt+1_pInt):(6_pInt*ns+3_pInt*nt), of)
|
||||
plasticState(ph)%state((5_pInt*ns+3_pInt*nt+nr+1_pInt):(6_pInt*ns+3_pInt*nt+nr), of)
|
||||
c = c + ns
|
||||
case (resolved_stress_slip_ID)
|
||||
j = 0_pInt
|
||||
|
@ -1822,7 +1830,7 @@ function constitutive_dislotwin_postResults(Tstar_v,Temperature,ipc,ip,el)
|
|||
c = c + ns
|
||||
case (threshold_stress_slip_ID)
|
||||
constitutive_dislotwin_postResults(c+1_pInt:c+ns) = &
|
||||
plasticState(ph)%state((6_pInt*ns+4_pInt*nt+1_pInt):(7_pInt*ns+4_pInt*nt), of)
|
||||
plasticState(ph)%state((6_pInt*ns+4_pInt*nt+nr+1_pInt):(7_pInt*ns+4_pInt*nt+nr), of)
|
||||
c = c + ns
|
||||
case (edge_dipole_distance_ID)
|
||||
j = 0_pInt
|
||||
|
@ -1834,9 +1842,9 @@ function constitutive_dislotwin_postResults(Tstar_v,Temperature,ipc,ip,el)
|
|||
(3.0_pReal*lattice_mu(ph)*constitutive_dislotwin_burgersPerSlipSystem(j,instance))/&
|
||||
(16.0_pReal*pi*abs(dot_product(Tstar_v,lattice_Sslip_v(:,1,index_myFamily+i,ph))))
|
||||
constitutive_dislotwin_postResults(c+j)=min(constitutive_dislotwin_postResults(c+j),&
|
||||
plasticState(ph)%state(5*ns+3*nt+j, of))
|
||||
plasticState(ph)%state(5*ns+3*nt+nr+j, of))
|
||||
! constitutive_dislotwin_postResults(c+j)=max(constitutive_dislotwin_postResults(c+j),&
|
||||
! plasticState(ph)%state(4*ns+2*nt+j, of))
|
||||
! plasticState(ph)%state(4*ns+2*nt+nr+j, of))
|
||||
enddo; enddo
|
||||
c = c + ns
|
||||
case (resolved_stress_shearband_ID)
|
||||
|
@ -1884,13 +1892,13 @@ function constitutive_dislotwin_postResults(Tstar_v,Temperature,ipc,ip,el)
|
|||
!* Resolved shear stress on slip system
|
||||
tau = dot_product(Tstar_v,lattice_Sslip_v(:,1,index_myFamily+i,ph))
|
||||
!* Stress ratios
|
||||
if((abs(tau)-plasticState(ph)%state(6*ns+4*nt+j, of)) > tol_math_check) then
|
||||
if((abs(tau)-plasticState(ph)%state(6*ns+4*nt+nr+j, of)) > tol_math_check) then
|
||||
!* Stress ratios
|
||||
StressRatio_p = ((abs(tau)-plasticState(ph)%state(6*ns+4*nt+j, of))/&
|
||||
StressRatio_p = ((abs(tau)-plasticState(ph)%state(6*ns+4*nt+nr+j, of))/&
|
||||
(constitutive_dislotwin_SolidSolutionStrength(instance)+&
|
||||
constitutive_dislotwin_tau_peierlsPerSlipFamily(f,instance)))&
|
||||
**constitutive_dislotwin_pPerSlipFamily(f,instance)
|
||||
StressRatio_pminus1 = ((abs(tau)-plasticState(ph)%state(6*ns+4*nt+j, of))/&
|
||||
StressRatio_pminus1 = ((abs(tau)-plasticState(ph)%state(6*ns+4*nt+nr+j, of))/&
|
||||
(constitutive_dislotwin_SolidSolutionStrength(instance)+&
|
||||
constitutive_dislotwin_tau_peierlsPerSlipFamily(f,instance)))&
|
||||
**(constitutive_dislotwin_pPerSlipFamily(f,instance)-1.0_pReal)
|
||||
|
@ -1918,7 +1926,7 @@ function constitutive_dislotwin_postResults(Tstar_v,Temperature,ipc,ip,el)
|
|||
!* Resolved shear stress on twin system
|
||||
tau = dot_product(Tstar_v,lattice_Stwin_v(:,index_myFamily+i,ph))
|
||||
!* Stress ratios
|
||||
StressRatio_r = (plasticState(ph)%state(7_pInt*ns+4_pInt*nt+j, of)/ &
|
||||
StressRatio_r = (plasticState(ph)%state(7_pInt*ns+4_pInt*nt+nr+j, of)/ &
|
||||
tau)**constitutive_dislotwin_rPerTwinFamily(f,instance)
|
||||
|
||||
!* Shear rates due to twin
|
||||
|
@ -1942,7 +1950,7 @@ function constitutive_dislotwin_postResults(Tstar_v,Temperature,ipc,ip,el)
|
|||
end select
|
||||
constitutive_dislotwin_postResults(c+j) = &
|
||||
(constitutive_dislotwin_MaxTwinFraction(instance)-sumf)*lattice_shearTwin(index_myFamily+i,ph)*&
|
||||
plasticState(ph)%state(7_pInt*ns+5_pInt*nt+j, of)*Ndot0*exp(-StressRatio_r)
|
||||
plasticState(ph)%state(7_pInt*ns+5_pInt*nt+nr+j, of)*Ndot0*exp(-StressRatio_r)
|
||||
endif
|
||||
|
||||
enddo ; enddo
|
||||
|
@ -1980,13 +1988,13 @@ function constitutive_dislotwin_postResults(Tstar_v,Temperature,ipc,ip,el)
|
|||
|
||||
!* Resolved shear stress on slip system
|
||||
tau = dot_product(Tstar_v,lattice_Sslip_v(:,1,index_myFamily+i,ph))
|
||||
if((abs(tau)-plasticState(ph)%state(6*ns+4*nt+j, of)) > tol_math_check) then
|
||||
if((abs(tau)-plasticState(ph)%state(6*ns+4*nt+nr+j, of)) > tol_math_check) then
|
||||
!* Stress ratios
|
||||
StressRatio_p = ((abs(tau)-plasticState(ph)%state(6*ns+4*nt+j, of))/&
|
||||
StressRatio_p = ((abs(tau)-plasticState(ph)%state(6*ns+4*nt+nr+j, of))/&
|
||||
(constitutive_dislotwin_SolidSolutionStrength(instance)+&
|
||||
constitutive_dislotwin_tau_peierlsPerSlipFamily(f,instance)))&
|
||||
**constitutive_dislotwin_pPerSlipFamily(f,instance)
|
||||
StressRatio_pminus1 = ((abs(tau)-plasticState(ph)%state(6*ns+4*nt+j, of))/&
|
||||
StressRatio_pminus1 = ((abs(tau)-plasticState(ph)%state(6*ns+4*nt+nr+j, of))/&
|
||||
(constitutive_dislotwin_SolidSolutionStrength(instance)+&
|
||||
constitutive_dislotwin_tau_peierlsPerSlipFamily(f,instance)))&
|
||||
**(constitutive_dislotwin_pPerSlipFamily(f,instance)-1.0_pReal)
|
||||
|
|
|
@ -1015,9 +1015,9 @@ subroutine lattice_init
|
|||
lattice_C66(6,6,section) = IO_floatValue(line,positions,2_pInt)
|
||||
case ('covera_ratio','c/a_ratio','c/a')
|
||||
CoverA(section) = IO_floatValue(line,positions,2_pInt)
|
||||
case ('aa', 'a_a', 'a_austenite')
|
||||
case ('aa', 'a_a', 'a_austenite', 'a_fcc')
|
||||
aA(section) = IO_floatValue(line,positions,2_pInt)
|
||||
case ('am', 'a_m', 'a_martensite')
|
||||
case ('am', 'a_m', 'a_martensite', 'a_bcc')
|
||||
aM(section) = IO_floatValue(line,positions,2_pInt)
|
||||
case ('cm', 'c_m', 'c_martensite')
|
||||
cM(section) = IO_floatValue(line,positions,2_pInt)
|
||||
|
|
Loading…
Reference in New Issue