missing rename: bcc -> cI
also, no CamelCase in disloTungsten + a few renames of variable names for more consistency
This commit is contained in:
parent
1d05edd7b3
commit
85d161f039
2
PRIVATE
2
PRIVATE
|
@ -1 +1 @@
|
||||||
Subproject commit 172ae90e9b8ae836071c7630e80fb0b70b044fed
|
Subproject commit 751f96155294e449b31c5d4913fea7bc02ce81b1
|
|
@ -30,10 +30,10 @@ submodule(constitutive) constitutive_mech
|
||||||
myPlasticity
|
myPlasticity
|
||||||
end function plastic_dislotwin_init
|
end function plastic_dislotwin_init
|
||||||
|
|
||||||
module function plastic_disloTungsten_init() result(myPlasticity)
|
module function plastic_dislotungsten_init() result(myPlasticity)
|
||||||
logical, dimension(:), allocatable :: &
|
logical, dimension(:), allocatable :: &
|
||||||
myPlasticity
|
myPlasticity
|
||||||
end function plastic_disloTungsten_init
|
end function plastic_dislotungsten_init
|
||||||
|
|
||||||
module function plastic_nonlocal_init() result(myPlasticity)
|
module function plastic_nonlocal_init() result(myPlasticity)
|
||||||
logical, dimension(:), allocatable :: &
|
logical, dimension(:), allocatable :: &
|
||||||
|
@ -94,7 +94,7 @@ submodule(constitutive) constitutive_mech
|
||||||
of
|
of
|
||||||
end subroutine plastic_dislotwin_LpAndItsTangent
|
end subroutine plastic_dislotwin_LpAndItsTangent
|
||||||
|
|
||||||
pure module subroutine plastic_disloTungsten_LpAndItsTangent(Lp,dLp_dMp,Mp,T,instance,of)
|
pure module subroutine plastic_dislotungsten_LpAndItsTangent(Lp,dLp_dMp,Mp,T,instance,of)
|
||||||
real(pReal), dimension(3,3), intent(out) :: &
|
real(pReal), dimension(3,3), intent(out) :: &
|
||||||
Lp !< plastic velocity gradient
|
Lp !< plastic velocity gradient
|
||||||
real(pReal), dimension(3,3,3,3), intent(out) :: &
|
real(pReal), dimension(3,3,3,3), intent(out) :: &
|
||||||
|
@ -107,7 +107,7 @@ submodule(constitutive) constitutive_mech
|
||||||
integer, intent(in) :: &
|
integer, intent(in) :: &
|
||||||
instance, &
|
instance, &
|
||||||
of
|
of
|
||||||
end subroutine plastic_disloTungsten_LpAndItsTangent
|
end subroutine plastic_dislotungsten_LpAndItsTangent
|
||||||
|
|
||||||
module subroutine plastic_nonlocal_LpAndItsTangent(Lp,dLp_dMp, &
|
module subroutine plastic_nonlocal_LpAndItsTangent(Lp,dLp_dMp, &
|
||||||
Mp,Temperature,instance,of,ip,el)
|
Mp,Temperature,instance,of,ip,el)
|
||||||
|
@ -136,11 +136,11 @@ submodule(constitutive) constitutive_mech
|
||||||
T
|
T
|
||||||
end subroutine plastic_dislotwin_dependentState
|
end subroutine plastic_dislotwin_dependentState
|
||||||
|
|
||||||
module subroutine plastic_disloTungsten_dependentState(instance,of)
|
module subroutine plastic_dislotungsten_dependentState(instance,of)
|
||||||
integer, intent(in) :: &
|
integer, intent(in) :: &
|
||||||
instance, &
|
instance, &
|
||||||
of
|
of
|
||||||
end subroutine plastic_disloTungsten_dependentState
|
end subroutine plastic_dislotungsten_dependentState
|
||||||
|
|
||||||
module subroutine plastic_nonlocal_dependentState(F, Fp, instance, of, ip, el)
|
module subroutine plastic_nonlocal_dependentState(F, Fp, instance, of, ip, el)
|
||||||
real(pReal), dimension(3,3), intent(in) :: &
|
real(pReal), dimension(3,3), intent(in) :: &
|
||||||
|
@ -173,10 +173,10 @@ submodule(constitutive) constitutive_mech
|
||||||
character(len=*), intent(in) :: group
|
character(len=*), intent(in) :: group
|
||||||
end subroutine plastic_dislotwin_results
|
end subroutine plastic_dislotwin_results
|
||||||
|
|
||||||
module subroutine plastic_disloTungsten_results(instance,group)
|
module subroutine plastic_dislotungsten_results(instance,group)
|
||||||
integer, intent(in) :: instance
|
integer, intent(in) :: instance
|
||||||
character(len=*), intent(in) :: group
|
character(len=*), intent(in) :: group
|
||||||
end subroutine plastic_disloTungsten_results
|
end subroutine plastic_dislotungsten_results
|
||||||
|
|
||||||
module subroutine plastic_nonlocal_results(instance,group)
|
module subroutine plastic_nonlocal_results(instance,group)
|
||||||
integer, intent(in) :: instance
|
integer, intent(in) :: instance
|
||||||
|
@ -255,7 +255,7 @@ module subroutine mech_init
|
||||||
where(plastic_phenopowerlaw_init()) phase_plasticity = PLASTICITY_PHENOPOWERLAW_ID
|
where(plastic_phenopowerlaw_init()) phase_plasticity = PLASTICITY_PHENOPOWERLAW_ID
|
||||||
where(plastic_kinehardening_init()) phase_plasticity = PLASTICITY_KINEHARDENING_ID
|
where(plastic_kinehardening_init()) phase_plasticity = PLASTICITY_KINEHARDENING_ID
|
||||||
where(plastic_dislotwin_init()) phase_plasticity = PLASTICITY_DISLOTWIN_ID
|
where(plastic_dislotwin_init()) phase_plasticity = PLASTICITY_DISLOTWIN_ID
|
||||||
where(plastic_disloTungsten_init()) phase_plasticity = PLASTICITY_DISLOTUNGSTEN_ID
|
where(plastic_dislotungsten_init()) phase_plasticity = PLASTICITY_DISLOTUNGSTEN_ID
|
||||||
where(plastic_nonlocal_init()) phase_plasticity = PLASTICITY_NONLOCAL_ID
|
where(plastic_nonlocal_init()) phase_plasticity = PLASTICITY_NONLOCAL_ID
|
||||||
|
|
||||||
do p = 1, phases%length
|
do p = 1, phases%length
|
||||||
|
@ -368,7 +368,7 @@ module subroutine constitutive_plastic_dependentState(F, Fp, ipc, ip, el)
|
||||||
case (PLASTICITY_DISLOTWIN_ID) plasticityType
|
case (PLASTICITY_DISLOTWIN_ID) plasticityType
|
||||||
call plastic_dislotwin_dependentState(temperature(ho)%p(tme),instance,of)
|
call plastic_dislotwin_dependentState(temperature(ho)%p(tme),instance,of)
|
||||||
case (PLASTICITY_DISLOTUNGSTEN_ID) plasticityType
|
case (PLASTICITY_DISLOTUNGSTEN_ID) plasticityType
|
||||||
call plastic_disloTungsten_dependentState(instance,of)
|
call plastic_dislotungsten_dependentState(instance,of)
|
||||||
case (PLASTICITY_NONLOCAL_ID) plasticityType
|
case (PLASTICITY_NONLOCAL_ID) plasticityType
|
||||||
call plastic_nonlocal_dependentState (F,Fp,instance,of,ip,el)
|
call plastic_nonlocal_dependentState (F,Fp,instance,of,ip,el)
|
||||||
end select plasticityType
|
end select plasticityType
|
||||||
|
@ -435,7 +435,7 @@ module subroutine constitutive_plastic_LpAndItsTangents(Lp, dLp_dS, dLp_dFi, &
|
||||||
call plastic_dislotwin_LpAndItsTangent(Lp,dLp_dMp,Mp,temperature(ho)%p(tme),instance,of)
|
call plastic_dislotwin_LpAndItsTangent(Lp,dLp_dMp,Mp,temperature(ho)%p(tme),instance,of)
|
||||||
|
|
||||||
case (PLASTICITY_DISLOTUNGSTEN_ID) plasticityType
|
case (PLASTICITY_DISLOTUNGSTEN_ID) plasticityType
|
||||||
call plastic_disloTungsten_LpAndItsTangent(Lp,dLp_dMp,Mp,temperature(ho)%p(tme),instance,of)
|
call plastic_dislotungsten_LpAndItsTangent(Lp,dLp_dMp,Mp,temperature(ho)%p(tme),instance,of)
|
||||||
|
|
||||||
end select plasticityType
|
end select plasticityType
|
||||||
|
|
||||||
|
@ -478,7 +478,7 @@ module subroutine plastic_results
|
||||||
call plastic_dislotwin_results(phase_plasticityInstance(p),group)
|
call plastic_dislotwin_results(phase_plasticityInstance(p),group)
|
||||||
|
|
||||||
case(PLASTICITY_DISLOTUNGSTEN_ID)
|
case(PLASTICITY_DISLOTUNGSTEN_ID)
|
||||||
call plastic_disloTungsten_results(phase_plasticityInstance(p),group)
|
call plastic_dislotungsten_results(phase_plasticityInstance(p),group)
|
||||||
|
|
||||||
case(PLASTICITY_NONLOCAL_ID)
|
case(PLASTICITY_NONLOCAL_ID)
|
||||||
call plastic_nonlocal_results(phase_plasticityInstance(p),group)
|
call plastic_nonlocal_results(phase_plasticityInstance(p),group)
|
||||||
|
|
|
@ -5,7 +5,7 @@
|
||||||
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
|
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
|
||||||
!> @brief crystal plasticity model for bcc metals, especially Tungsten
|
!> @brief crystal plasticity model for bcc metals, especially Tungsten
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
submodule(constitutive:constitutive_mech) plastic_disloTungsten
|
submodule(constitutive:constitutive_mech) plastic_dislotungsten
|
||||||
|
|
||||||
real(pReal), parameter :: &
|
real(pReal), parameter :: &
|
||||||
kB = 1.38e-23_pReal !< Boltzmann constant in J/Kelvin
|
kB = 1.38e-23_pReal !< Boltzmann constant in J/Kelvin
|
||||||
|
@ -74,7 +74,7 @@ contains
|
||||||
!> @brief Perform module initialization.
|
!> @brief Perform module initialization.
|
||||||
!> @details reads in material parameters, allocates arrays, and does sanity checks
|
!> @details reads in material parameters, allocates arrays, and does sanity checks
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
module function plastic_disloTungsten_init() result(myPlasticity)
|
module function plastic_dislotungsten_init() result(myPlasticity)
|
||||||
|
|
||||||
logical, dimension(:), allocatable :: myPlasticity
|
logical, dimension(:), allocatable :: myPlasticity
|
||||||
integer :: &
|
integer :: &
|
||||||
|
@ -99,7 +99,7 @@ module function plastic_disloTungsten_init() result(myPlasticity)
|
||||||
|
|
||||||
print'(/,a)', ' <<<+- plastic_dislotungsten init -+>>>'
|
print'(/,a)', ' <<<+- plastic_dislotungsten init -+>>>'
|
||||||
|
|
||||||
myPlasticity = plastic_active('disloTungsten')
|
myPlasticity = plastic_active('dislotungsten')
|
||||||
Ninstances = count(myPlasticity)
|
Ninstances = count(myPlasticity)
|
||||||
print'(a,i2)', ' # instances: ',Ninstances; flush(IO_STDOUT)
|
print'(a,i2)', ' # instances: ',Ninstances; flush(IO_STDOUT)
|
||||||
if(Ninstances == 0) return
|
if(Ninstances == 0) return
|
||||||
|
@ -142,7 +142,7 @@ module function plastic_disloTungsten_init() result(myPlasticity)
|
||||||
prm%P_sl = lattice_SchmidMatrix_slip(N_sl,phase%get_asString('lattice'),&
|
prm%P_sl = lattice_SchmidMatrix_slip(N_sl,phase%get_asString('lattice'),&
|
||||||
phase%get_asFloat('c/a',defaultVal=0.0_pReal))
|
phase%get_asFloat('c/a',defaultVal=0.0_pReal))
|
||||||
|
|
||||||
if(trim(phase%get_asString('lattice')) == 'bcc') then
|
if(trim(phase%get_asString('lattice')) == 'cI') then
|
||||||
a = pl%get_asFloats('a_nonSchmid',defaultVal = emptyRealArray)
|
a = pl%get_asFloats('a_nonSchmid',defaultVal = emptyRealArray)
|
||||||
prm%nonSchmid_pos = lattice_nonSchmidMatrix(N_sl,a,+1)
|
prm%nonSchmid_pos = lattice_nonSchmidMatrix(N_sl,a,+1)
|
||||||
prm%nonSchmid_neg = lattice_nonSchmidMatrix(N_sl,a,-1)
|
prm%nonSchmid_neg = lattice_nonSchmidMatrix(N_sl,a,-1)
|
||||||
|
@ -180,7 +180,7 @@ module function plastic_disloTungsten_init() result(myPlasticity)
|
||||||
prm%f_at = pl%get_asFloat('f_at') * prm%b_sl**3.0_pReal
|
prm%f_at = pl%get_asFloat('f_at') * prm%b_sl**3.0_pReal
|
||||||
prm%D_a = pl%get_asFloat('D_a') * prm%b_sl
|
prm%D_a = pl%get_asFloat('D_a') * prm%b_sl
|
||||||
|
|
||||||
prm%dipoleformation = pl%get_asBool('dipole_formation_factor', defaultVal = .true.)
|
prm%dipoleformation = .not. pl%get_asBool('no_dipole_formation', defaultVal = .false.)
|
||||||
|
|
||||||
! expand: family => system
|
! expand: family => system
|
||||||
rho_mob_0 = math_expand(rho_mob_0, N_sl)
|
rho_mob_0 = math_expand(rho_mob_0, N_sl)
|
||||||
|
@ -262,17 +262,17 @@ module function plastic_disloTungsten_init() result(myPlasticity)
|
||||||
|
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
! exit if any parameter is out of range
|
! exit if any parameter is out of range
|
||||||
if (extmsg /= '') call IO_error(211,ext_msg=trim(extmsg)//'(disloTungsten)')
|
if (extmsg /= '') call IO_error(211,ext_msg=trim(extmsg)//'(dislotungsten)')
|
||||||
|
|
||||||
enddo
|
enddo
|
||||||
|
|
||||||
end function plastic_disloTungsten_init
|
end function plastic_dislotungsten_init
|
||||||
|
|
||||||
|
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
!> @brief Calculate plastic velocity gradient and its tangent.
|
!> @brief Calculate plastic velocity gradient and its tangent.
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
pure module subroutine plastic_disloTungsten_LpAndItsTangent(Lp,dLp_dMp, &
|
pure module subroutine plastic_dislotungsten_LpAndItsTangent(Lp,dLp_dMp, &
|
||||||
Mp,T,instance,of)
|
Mp,T,instance,of)
|
||||||
real(pReal), dimension(3,3), intent(out) :: &
|
real(pReal), dimension(3,3), intent(out) :: &
|
||||||
Lp !< plastic velocity gradient
|
Lp !< plastic velocity gradient
|
||||||
|
@ -309,13 +309,13 @@ pure module subroutine plastic_disloTungsten_LpAndItsTangent(Lp,dLp_dMp, &
|
||||||
|
|
||||||
end associate
|
end associate
|
||||||
|
|
||||||
end subroutine plastic_disloTungsten_LpAndItsTangent
|
end subroutine plastic_dislotungsten_LpAndItsTangent
|
||||||
|
|
||||||
|
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
!> @brief Calculate the rate of change of microstructure.
|
!> @brief Calculate the rate of change of microstructure.
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
module subroutine plastic_disloTungsten_dotState(Mp,T,instance,of)
|
module subroutine plastic_dislotungsten_dotState(Mp,T,instance,of)
|
||||||
|
|
||||||
real(pReal), dimension(3,3), intent(in) :: &
|
real(pReal), dimension(3,3), intent(in) :: &
|
||||||
Mp !< Mandel stress
|
Mp !< Mandel stress
|
||||||
|
@ -369,13 +369,13 @@ module subroutine plastic_disloTungsten_dotState(Mp,T,instance,of)
|
||||||
|
|
||||||
end associate
|
end associate
|
||||||
|
|
||||||
end subroutine plastic_disloTungsten_dotState
|
end subroutine plastic_dislotungsten_dotState
|
||||||
|
|
||||||
|
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
!> @brief Calculate derived quantities from state.
|
!> @brief Calculate derived quantities from state.
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
module subroutine plastic_disloTungsten_dependentState(instance,of)
|
module subroutine plastic_dislotungsten_dependentState(instance,of)
|
||||||
|
|
||||||
integer, intent(in) :: &
|
integer, intent(in) :: &
|
||||||
instance, &
|
instance, &
|
||||||
|
@ -394,13 +394,13 @@ module subroutine plastic_disloTungsten_dependentState(instance,of)
|
||||||
|
|
||||||
end associate
|
end associate
|
||||||
|
|
||||||
end subroutine plastic_disloTungsten_dependentState
|
end subroutine plastic_dislotungsten_dependentState
|
||||||
|
|
||||||
|
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
!> @brief Write results to HDF5 output file.
|
!> @brief Write results to HDF5 output file.
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
module subroutine plastic_disloTungsten_results(instance,group)
|
module subroutine plastic_dislotungsten_results(instance,group)
|
||||||
|
|
||||||
integer, intent(in) :: instance
|
integer, intent(in) :: instance
|
||||||
character(len=*), intent(in) :: group
|
character(len=*), intent(in) :: group
|
||||||
|
@ -429,7 +429,7 @@ module subroutine plastic_disloTungsten_results(instance,group)
|
||||||
enddo outputsLoop
|
enddo outputsLoop
|
||||||
end associate
|
end associate
|
||||||
|
|
||||||
end subroutine plastic_disloTungsten_results
|
end subroutine plastic_dislotungsten_results
|
||||||
|
|
||||||
|
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
|
@ -547,4 +547,4 @@ pure subroutine kinetics(Mp,T,instance,of, &
|
||||||
|
|
||||||
end subroutine kinetics
|
end subroutine kinetics
|
||||||
|
|
||||||
end submodule plastic_disloTungsten
|
end submodule plastic_dislotungsten
|
||||||
|
|
|
@ -301,7 +301,7 @@ module function plastic_dislotwin_init() result(myPlasticity)
|
||||||
prm%r = math_expand(prm%r,N_tw)
|
prm%r = math_expand(prm%r,N_tw)
|
||||||
|
|
||||||
! sanity checks
|
! sanity checks
|
||||||
if ( prm%x_c_tw < 0.0_pReal) extmsg = trim(extmsg)//' x_c_twin'
|
if ( prm%x_c_tw < 0.0_pReal) extmsg = trim(extmsg)//' x_c_tw'
|
||||||
if ( prm%L_tw < 0.0_pReal) extmsg = trim(extmsg)//' L_tw'
|
if ( prm%L_tw < 0.0_pReal) extmsg = trim(extmsg)//' L_tw'
|
||||||
if ( prm%i_tw < 0.0_pReal) extmsg = trim(extmsg)//' i_tw'
|
if ( prm%i_tw < 0.0_pReal) extmsg = trim(extmsg)//' i_tw'
|
||||||
if (any(prm%b_tw < 0.0_pReal)) extmsg = trim(extmsg)//' b_tw'
|
if (any(prm%b_tw < 0.0_pReal)) extmsg = trim(extmsg)//' b_tw'
|
||||||
|
@ -332,15 +332,15 @@ module function plastic_dislotwin_init() result(myPlasticity)
|
||||||
prm%h_tr_tr = lattice_interaction_TransByTrans(N_tr,pl%get_asFloats('h_tr_tr'), &
|
prm%h_tr_tr = lattice_interaction_TransByTrans(N_tr,pl%get_asFloats('h_tr_tr'), &
|
||||||
phase%get_asString('lattice'))
|
phase%get_asString('lattice'))
|
||||||
|
|
||||||
prm%C66_tr = lattice_C66_trans(N_tr,prm%C66,pl%get_asString('trans_lattice_structure'), &
|
prm%C66_tr = lattice_C66_trans(N_tr,prm%C66,pl%get_asString('lattice_tr'), &
|
||||||
0.0_pReal, &
|
0.0_pReal, &
|
||||||
pl%get_asFloat('a_bcc', defaultVal=0.0_pReal), &
|
pl%get_asFloat('a_cI', defaultVal=0.0_pReal), &
|
||||||
pl%get_asFloat('a_fcc', defaultVal=0.0_pReal))
|
pl%get_asFloat('a_cF', defaultVal=0.0_pReal))
|
||||||
|
|
||||||
prm%P_tr = lattice_SchmidMatrix_trans(N_tr,pl%get_asString('trans_lattice_structure'), &
|
prm%P_tr = lattice_SchmidMatrix_trans(N_tr,pl%get_asString('lattice_tr'), &
|
||||||
0.0_pReal, &
|
0.0_pReal, &
|
||||||
pl%get_asFloat('a_bcc', defaultVal=0.0_pReal), &
|
pl%get_asFloat('a_cI', defaultVal=0.0_pReal), &
|
||||||
pl%get_asFloat('a_fcc', defaultVal=0.0_pReal))
|
pl%get_asFloat('a_cF', defaultVal=0.0_pReal))
|
||||||
|
|
||||||
if (lattice_structure(p) /= lattice_FCC_ID) then
|
if (lattice_structure(p) /= lattice_FCC_ID) then
|
||||||
prm%dot_N_0_tr = pl%get_asFloats('dot_N_0_tr')
|
prm%dot_N_0_tr = pl%get_asFloats('dot_N_0_tr')
|
||||||
|
@ -352,7 +352,7 @@ module function plastic_dislotwin_init() result(myPlasticity)
|
||||||
prm%s = math_expand(prm%s,N_tr)
|
prm%s = math_expand(prm%s,N_tr)
|
||||||
|
|
||||||
! sanity checks
|
! sanity checks
|
||||||
if ( prm%x_c_tr < 0.0_pReal) extmsg = trim(extmsg)//' x_c_trans'
|
if ( prm%x_c_tr < 0.0_pReal) extmsg = trim(extmsg)//' x_c_tr'
|
||||||
if ( prm%L_tr < 0.0_pReal) extmsg = trim(extmsg)//' L_tr'
|
if ( prm%L_tr < 0.0_pReal) extmsg = trim(extmsg)//' L_tr'
|
||||||
if ( prm%i_tr < 0.0_pReal) extmsg = trim(extmsg)//' i_tr'
|
if ( prm%i_tr < 0.0_pReal) extmsg = trim(extmsg)//' i_tr'
|
||||||
if (any(prm%t_tr < 0.0_pReal)) extmsg = trim(extmsg)//' t_tr'
|
if (any(prm%t_tr < 0.0_pReal)) extmsg = trim(extmsg)//' t_tr'
|
||||||
|
|
|
@ -125,7 +125,7 @@ module function plastic_kinehardening_init() result(myPlasticity)
|
||||||
prm%P = lattice_SchmidMatrix_slip(N_sl,phase%get_asString('lattice'),&
|
prm%P = lattice_SchmidMatrix_slip(N_sl,phase%get_asString('lattice'),&
|
||||||
phase%get_asFloat('c/a',defaultVal=0.0_pReal))
|
phase%get_asFloat('c/a',defaultVal=0.0_pReal))
|
||||||
|
|
||||||
if(trim(phase%get_asString('lattice')) == 'bcc') then
|
if(trim(phase%get_asString('lattice')) == 'cI') then
|
||||||
a = pl%get_asFloats('a_nonSchmid',defaultVal = emptyRealArray)
|
a = pl%get_asFloats('a_nonSchmid',defaultVal = emptyRealArray)
|
||||||
if(size(a) > 0) prm%nonSchmidActive = .true.
|
if(size(a) > 0) prm%nonSchmidActive = .true.
|
||||||
prm%nonSchmid_pos = lattice_nonSchmidMatrix(N_sl,a,+1)
|
prm%nonSchmid_pos = lattice_nonSchmidMatrix(N_sl,a,+1)
|
||||||
|
|
|
@ -244,7 +244,7 @@ module function plastic_nonlocal_init() result(myPlasticity)
|
||||||
prm%Schmid = lattice_SchmidMatrix_slip(ini%N_sl,phase%get_asString('lattice'),&
|
prm%Schmid = lattice_SchmidMatrix_slip(ini%N_sl,phase%get_asString('lattice'),&
|
||||||
phase%get_asFloat('c/a',defaultVal=0.0_pReal))
|
phase%get_asFloat('c/a',defaultVal=0.0_pReal))
|
||||||
|
|
||||||
if(trim(phase%get_asString('lattice')) == 'bcc') then
|
if(trim(phase%get_asString('lattice')) == 'cI') then
|
||||||
a = pl%get_asFloats('a_nonSchmid',defaultVal = emptyRealArray)
|
a = pl%get_asFloats('a_nonSchmid',defaultVal = emptyRealArray)
|
||||||
if(size(a) > 0) prm%nonSchmidActive = .true.
|
if(size(a) > 0) prm%nonSchmidActive = .true.
|
||||||
prm%nonSchmid_pos = lattice_nonSchmidMatrix(ini%N_sl,a,+1)
|
prm%nonSchmid_pos = lattice_nonSchmidMatrix(ini%N_sl,a,+1)
|
||||||
|
@ -362,7 +362,7 @@ module function plastic_nonlocal_init() result(myPlasticity)
|
||||||
if (prm%nu_a <= 0.0_pReal) extmsg = trim(extmsg)//' nu_a'
|
if (prm%nu_a <= 0.0_pReal) extmsg = trim(extmsg)//' nu_a'
|
||||||
if (prm%w <= 0.0_pReal) extmsg = trim(extmsg)//' w'
|
if (prm%w <= 0.0_pReal) extmsg = trim(extmsg)//' w'
|
||||||
if (prm%D_0 < 0.0_pReal) extmsg = trim(extmsg)//' D_0'
|
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%V_at <= 0.0_pReal) extmsg = trim(extmsg)//' V_at' ! ToDo: in dislotungsten, the atomic volume is given as a factor
|
||||||
|
|
||||||
if (prm%rho_min < 0.0_pReal) extmsg = trim(extmsg)//' rho_min'
|
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%rho_significant < 0.0_pReal) extmsg = trim(extmsg)//' rho_significant'
|
||||||
|
|
|
@ -120,7 +120,7 @@ module function plastic_phenopowerlaw_init() result(myPlasticity)
|
||||||
prm%P_sl = lattice_SchmidMatrix_slip(N_sl,phase%get_asString('lattice'),&
|
prm%P_sl = lattice_SchmidMatrix_slip(N_sl,phase%get_asString('lattice'),&
|
||||||
phase%get_asFloat('c/a',defaultVal=0.0_pReal))
|
phase%get_asFloat('c/a',defaultVal=0.0_pReal))
|
||||||
|
|
||||||
if(phase%get_asString('lattice') == 'bcc') then
|
if(phase%get_asString('lattice') == 'cI') then
|
||||||
a = pl%get_asFloats('a_nonSchmid',defaultVal=emptyRealArray)
|
a = pl%get_asFloats('a_nonSchmid',defaultVal=emptyRealArray)
|
||||||
if(size(a) > 0) prm%nonSchmidActive = .true.
|
if(size(a) > 0) prm%nonSchmidActive = .true.
|
||||||
prm%nonSchmid_pos = lattice_nonSchmidMatrix(N_sl,a,+1)
|
prm%nonSchmid_pos = lattice_nonSchmidMatrix(N_sl,a,+1)
|
||||||
|
|
Loading…
Reference in New Issue