improved parameter names

This commit is contained in:
Sharan Roongta 2020-09-23 00:55:19 +02:00
parent 881167a37c
commit de3e13df88
5 changed files with 29 additions and 29 deletions

@ -1 +1 @@
Subproject commit 555f3e01f2b5cf43ade1bd48423b890adca21771
Subproject commit f7b371f12950398b61b1ab41fe126c6045e4e686

View File

@ -21,7 +21,7 @@ submodule(constitutive:constitutive_plastic) plastic_disloTungsten
D_a, &
i_sl, & !< Adj. parameter for distance between 2 forest dislocations
f_at, & !< factor to calculate atomic volume
tau_peierls, & !< Peierls stress
tau_Peierls, & !< Peierls stress
!* mobility law parameters
Q_s, & !< activation energy for glide [J]
v_0, & !< dislocation velocity prefactor [m/s]
@ -142,7 +142,7 @@ module function plastic_disloTungsten_init() result(myPlasticity)
phase%get_asFloat('c/a',defaultVal=0.0_pReal))
if(trim(phase%get_asString('lattice')) == 'bcc') then
a = pl%get_asFloats('nonSchmid_coefficients',defaultVal = emptyRealArray)
a = pl%get_asFloats('a_nonSchmid',defaultVal = emptyRealArray)
prm%nonSchmid_pos = lattice_nonSchmidMatrix(N_sl,a,+1)
prm%nonSchmid_neg = lattice_nonSchmidMatrix(N_sl,a,-1)
else
@ -163,7 +163,7 @@ module function plastic_disloTungsten_init() result(myPlasticity)
prm%Q_s = pl%get_asFloats('Q_s', requiredSize=size(N_sl))
prm%i_sl = pl%get_asFloats('i_sl', requiredSize=size(N_sl))
prm%tau_peierls = pl%get_asFloats('tau_peierls', requiredSize=size(N_sl))
prm%tau_Peierls = pl%get_asFloats('tau_Peierls', requiredSize=size(N_sl))
prm%p = pl%get_asFloats('p_sl', requiredSize=size(N_sl), &
defaultVal=[(1.0_pReal,i=1,size(N_sl))])
prm%q = pl%get_asFloats('q_sl', requiredSize=size(N_sl), &
@ -191,7 +191,7 @@ module function plastic_disloTungsten_init() result(myPlasticity)
prm%h = math_expand(prm%h, N_sl)
prm%w = math_expand(prm%w, N_sl)
prm%omega = math_expand(prm%omega, N_sl)
prm%tau_peierls = math_expand(prm%tau_peierls, N_sl)
prm%tau_Peierls = math_expand(prm%tau_Peierls, N_sl)
prm%v_0 = math_expand(prm%v_0, N_sl)
prm%B = math_expand(prm%B, N_sl)
prm%i_sl = math_expand(prm%i_sl, N_sl)
@ -206,13 +206,13 @@ module function plastic_disloTungsten_init() result(myPlasticity)
if (any(prm%v_0 < 0.0_pReal)) extmsg = trim(extmsg)//' v_0'
if (any(prm%b_sl <= 0.0_pReal)) extmsg = trim(extmsg)//' b_sl'
if (any(prm%Q_s <= 0.0_pReal)) extmsg = trim(extmsg)//' Q_s'
if (any(prm%tau_peierls < 0.0_pReal)) extmsg = trim(extmsg)//' tau_peierls'
if (any(prm%tau_Peierls < 0.0_pReal)) extmsg = trim(extmsg)//' tau_Peierls'
if (any(prm%D_a <= 0.0_pReal)) extmsg = trim(extmsg)//' D_a or b_sl'
if (any(prm%f_at <= 0.0_pReal)) extmsg = trim(extmsg)//' f_at or b_sl'
else slipActive
rho_mob_0= emptyRealArray; rho_dip_0 = emptyRealArray
allocate(prm%b_sl,prm%D_a,prm%i_sl,prm%f_at,prm%tau_peierls, &
allocate(prm%b_sl,prm%D_a,prm%i_sl,prm%f_at,prm%tau_Peierls, &
prm%Q_s,prm%v_0,prm%p,prm%q,prm%B,prm%h,prm%w,prm%omega, &
source = emptyRealArray)
allocate(prm%forestProjection(0,0))
@ -482,7 +482,7 @@ pure subroutine kinetics(Mp,T,instance,of, &
effectiveLength => dst%Lambda_sl(:,of) - prm%w)
significantPositiveTau: where(abs(tau_pos)-dst%threshold_stress(:,of) > tol_math_check)
StressRatio = (abs(tau_pos)-dst%threshold_stress(:,of))/prm%tau_peierls
StressRatio = (abs(tau_pos)-dst%threshold_stress(:,of))/prm%tau_Peierls
StressRatio_p = StressRatio** prm%p
StressRatio_pminus1 = StressRatio**(prm%p-1.0_pReal)
needsGoodName = exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q)
@ -500,7 +500,7 @@ pure subroutine kinetics(Mp,T,instance,of, &
if (present(ddot_gamma_dtau_pos)) then
significantPositiveTau2: where(abs(tau_pos)-dst%threshold_stress(:,of) > tol_math_check)
dtn = -1.0_pReal * t_n * BoltzmannRatio * prm%p * prm%q * (1.0_pReal-StressRatio_p)**(prm%q - 1.0_pReal) &
* (StressRatio)**(prm%p - 1.0_pReal) / prm%tau_peierls
* (StressRatio)**(prm%p - 1.0_pReal) / prm%tau_Peierls
dtk = -1.0_pReal * t_k / tau_pos
dvel = -1.0_pReal * prm%h * (dtk + dtn) / (t_n + t_k)**2.0_pReal
@ -512,7 +512,7 @@ pure subroutine kinetics(Mp,T,instance,of, &
endif
significantNegativeTau: where(abs(tau_neg)-dst%threshold_stress(:,of) > tol_math_check)
StressRatio = (abs(tau_neg)-dst%threshold_stress(:,of))/prm%tau_peierls
StressRatio = (abs(tau_neg)-dst%threshold_stress(:,of))/prm%tau_Peierls
StressRatio_p = StressRatio** prm%p
StressRatio_pminus1 = StressRatio**(prm%p-1.0_pReal)
needsGoodName = exp(-BoltzmannRatio*(1-StressRatio_p) ** prm%q)
@ -530,7 +530,7 @@ pure subroutine kinetics(Mp,T,instance,of, &
if (present(ddot_gamma_dtau_neg)) then
significantNegativeTau2: where(abs(tau_neg)-dst%threshold_stress(:,of) > tol_math_check)
dtn = -1.0_pReal * t_n * BoltzmannRatio * prm%p * prm%q * (1.0_pReal-StressRatio_p)**(prm%q - 1.0_pReal) &
* (StressRatio)**(prm%p - 1.0_pReal) / prm%tau_peierls
* (StressRatio)**(prm%p - 1.0_pReal) / prm%tau_Peierls
dtk = -1.0_pReal * t_k / tau_neg
dvel = -1.0_pReal * prm%h * (dtk + dtn) / (t_n + t_k)**2.0_pReal

View File

@ -125,7 +125,7 @@ module function plastic_kinehardening_init() result(myPlasticity)
phase%get_asFloat('c/a',defaultVal=0.0_pReal))
if(trim(phase%get_asString('lattice')) == 'bcc') then
a = pl%get_asFloats('nonSchmid_coefficients',defaultVal = emptyRealArray)
a = pl%get_asFloats('a_nonSchmid',defaultVal = emptyRealArray)
if(size(a) > 0) prm%nonSchmidActive = .true.
prm%nonSchmid_pos = lattice_nonSchmidMatrix(N_sl,a,+1)
prm%nonSchmid_neg = lattice_nonSchmidMatrix(N_sl,a,-1)

View File

@ -67,8 +67,8 @@ submodule(constitutive:constitutive_plastic) plastic_nonlocal
Q_cl, & !< activation enthalpy for diffusion
atol_rho, & !< absolute tolerance for dislocation density in state integration
rho_significant, & !< density considered significant
rho_num_significant, & !< number of dislocations considered significant
w, & !< width of a doubkle kink in multiples of the burgers vector length b
rho_min, & !< number of dislocations considered significant
w, & !< width of a doubkle kink in multiples of the Burgers vector length b
Q_sol, & !< activation energy for solid solution in J
f_sol, & !< solid solution obstacle size in multiples of the burgers vector length
c_sol, & !< concentration of solid solution in atomic parts
@ -87,8 +87,8 @@ submodule(constitutive:constitutive_plastic) plastic_nonlocal
real(pReal), dimension(:), allocatable :: &
d_ed, & !< minimum stable edge dipole height
d_sc, & !< minimum stable screw dipole height
tau_peierls_ed, &
tau_peierls_sc, &
tau_Peierls_ed, &
tau_Peierls_sc, &
i_sl, & !< mean free path prefactor for each
b_sl !< absolute length of burgers vector [m]
real(pReal), dimension(:,:), allocatable :: &
@ -244,7 +244,7 @@ module function plastic_nonlocal_init() result(myPlasticity)
phase%get_asFloat('c/a',defaultVal=0.0_pReal))
if(trim(phase%get_asString('lattice')) == 'bcc') then
a = pl%get_asFloats('nonSchmid_coefficients',defaultVal = emptyRealArray)
a = pl%get_asFloats('a_nonSchmid',defaultVal = emptyRealArray)
if(size(a) > 0) prm%nonSchmidActive = .true.
prm%nonSchmid_pos = lattice_nonSchmidMatrix(ini%N_sl,a,+1)
prm%nonSchmid_neg = lattice_nonSchmidMatrix(ini%N_sl,a,-1)
@ -300,16 +300,16 @@ module function plastic_nonlocal_init() result(myPlasticity)
prm%minDipoleHeight(:,1) = prm%d_ed
prm%minDipoleHeight(:,2) = prm%d_sc
prm%tau_peierls_ed = pl%get_asFloats('tau_peierls_ed', requiredSize=size(ini%N_sl))
prm%tau_peierls_sc = pl%get_asFloats('tau_peierls_sc', requiredSize=size(ini%N_sl))
prm%tau_peierls_ed = math_expand(prm%tau_peierls_ed,ini%N_sl)
prm%tau_peierls_sc = math_expand(prm%tau_peierls_sc,ini%N_sl)
prm%tau_Peierls_ed = pl%get_asFloats('tau_Peierls_ed', requiredSize=size(ini%N_sl))
prm%tau_Peierls_sc = pl%get_asFloats('tau_Peierls_sc', requiredSize=size(ini%N_sl))
prm%tau_Peierls_ed = math_expand(prm%tau_Peierls_ed,ini%N_sl)
prm%tau_Peierls_sc = math_expand(prm%tau_Peierls_sc,ini%N_sl)
allocate(prm%peierlsstress(prm%sum_N_sl,2))
prm%peierlsstress(:,1) = prm%tau_peierls_ed
prm%peierlsstress(:,2) = prm%tau_peierls_sc
prm%peierlsstress(:,1) = prm%tau_Peierls_ed
prm%peierlsstress(:,2) = prm%tau_Peierls_sc
prm%rho_significant = pl%get_asFloat('rho_significant')
prm%rho_num_significant = pl%get_asFloat('rho_num_significant', 0.0_pReal)
prm%rho_min = pl%get_asFloat('rho_min', 0.0_pReal)
prm%f_c = pl%get_asFloat('f_c',defaultVal=2.0_pReal)
prm%V_at = pl%get_asFloat('V_at')
@ -363,7 +363,7 @@ module function plastic_nonlocal_init() result(myPlasticity)
if (prm%D_0 < 0.0_pReal) extmsg = trim(extmsg)//' D_0'
if (prm%V_at <= 0.0_pReal) extmsg = trim(extmsg)//' V_at' ! ToDo: in disloTungsten, the atomic volume is given as a factor
if (prm%rho_num_significant < 0.0_pReal) extmsg = trim(extmsg)//' rho_num_significant'
if (prm%rho_min < 0.0_pReal) extmsg = trim(extmsg)//' rho_min'
if (prm%rho_significant < 0.0_pReal) extmsg = trim(extmsg)//' rho_significant'
if (prm%atol_rho < 0.0_pReal) extmsg = trim(extmsg)//' atol_rho'
if (prm%f_c < 0.0_pReal) extmsg = trim(extmsg)//' f_c'
@ -1334,7 +1334,7 @@ function rhoDotFlux(F,Fp,timestep, instance,of,ip,el)
neighbor_rhoSgl0(s,t) = max(plasticState(np)%state0(iRhoU(s,t,neighbor_instance),no),0.0_pReal)
endforall
where (neighbor_rhoSgl0 * IPvolume(neighbor_ip,neighbor_el) ** 0.667_pReal < prm%rho_num_significant &
where (neighbor_rhoSgl0 * IPvolume(neighbor_ip,neighbor_el) ** 0.667_pReal < prm%rho_min &
.or. neighbor_rhoSgl0 < prm%rho_significant) &
neighbor_rhoSgl0 = 0.0_pReal
normal_neighbor2me_defConf = math_det33(Favg) * matmul(math_inv33(transpose(Favg)), &
@ -1805,7 +1805,7 @@ pure function getRho(instance,of,ip,el)
getRho(:,mob) = max(getRho(:,mob),0.0_pReal)
getRho(:,dip) = max(getRho(:,dip),0.0_pReal)
where(abs(getRho) < max(prm%rho_num_significant/IPvolume(ip,el)**(2.0_pReal/3.0_pReal),prm%rho_significant)) &
where(abs(getRho) < max(prm%rho_min/IPvolume(ip,el)**(2.0_pReal/3.0_pReal),prm%rho_significant)) &
getRho = 0.0_pReal
end associate
@ -1830,7 +1830,7 @@ pure function getRho0(instance,of,ip,el)
getRho0(:,mob) = max(getRho0(:,mob),0.0_pReal)
getRho0(:,dip) = max(getRho0(:,dip),0.0_pReal)
where(abs(getRho0) < max(prm%rho_num_significant/IPvolume(ip,el)**(2.0_pReal/3.0_pReal),prm%rho_significant)) &
where(abs(getRho0) < max(prm%rho_min/IPvolume(ip,el)**(2.0_pReal/3.0_pReal),prm%rho_significant)) &
getRho0 = 0.0_pReal
end associate

View File

@ -120,7 +120,7 @@ module function plastic_phenopowerlaw_init() result(myPlasticity)
phase%get_asFloat('c/a',defaultVal=0.0_pReal))
if(phase%get_asString('lattice') == 'bcc') then
a = pl%get_asFloats('nonSchmid_coefficients',defaultVal=emptyRealArray)
a = pl%get_asFloats('a_nonSchmid',defaultVal=emptyRealArray)
if(size(a) > 0) prm%nonSchmidActive = .true.
prm%nonSchmid_pos = lattice_nonSchmidMatrix(N_sl,a,+1)
prm%nonSchmid_neg = lattice_nonSchmidMatrix(N_sl,a,-1)