put private functions at the end

for easy separation
This commit is contained in:
Martin Diehl 2018-12-21 14:56:32 +01:00
parent fc171f388a
commit dcd22ccb6a
1 changed files with 103 additions and 99 deletions

View File

@ -420,9 +420,9 @@ subroutine plastic_phenopowerlaw_LpAndItsTangent(Lp,dLp_dMp,Mp,instance,of)
Lp = 0.0_pReal Lp = 0.0_pReal
dLp_dMp = 0.0_pReal dLp_dMp = 0.0_pReal
associate(prm => param(instance), stt => state(instance)) associate(prm => param(instance))
call kinetics_slip(prm,stt,of,Mp,gdot_slip_pos,gdot_slip_neg,dgdot_dtauslip_pos,dgdot_dtauslip_neg) call kinetics_slip(Mp,instance,of,gdot_slip_pos,gdot_slip_neg,dgdot_dtauslip_pos,dgdot_dtauslip_neg)
slipSystems: do i = 1_pInt, prm%totalNslip slipSystems: do i = 1_pInt, prm%totalNslip
Lp = Lp + (gdot_slip_pos(i)+gdot_slip_neg(i))*prm%Schmid_slip(1:3,1:3,i) Lp = Lp + (gdot_slip_pos(i)+gdot_slip_neg(i))*prm%Schmid_slip(1:3,1:3,i)
forall (k=1_pInt:3_pInt,l=1_pInt:3_pInt,m=1_pInt:3_pInt,n=1_pInt:3_pInt) & forall (k=1_pInt:3_pInt,l=1_pInt:3_pInt,m=1_pInt:3_pInt,n=1_pInt:3_pInt) &
@ -431,7 +431,7 @@ subroutine plastic_phenopowerlaw_LpAndItsTangent(Lp,dLp_dMp,Mp,instance,of)
+ dgdot_dtauslip_neg(i) * prm%Schmid_slip(k,l,i) * prm%nonSchmid_neg(m,n,i) + dgdot_dtauslip_neg(i) * prm%Schmid_slip(k,l,i) * prm%nonSchmid_neg(m,n,i)
enddo slipSystems enddo slipSystems
call kinetics_twin(prm,stt,of,Mp,gdot_twin,dgdot_dtautwin) call kinetics_twin(Mp,instance,of,gdot_twin,dgdot_dtautwin)
twinSystems: do i = 1_pInt, prm%totalNtwin twinSystems: do i = 1_pInt, prm%totalNtwin
Lp = Lp + gdot_twin(i)*prm%Schmid_twin(1:3,1:3,i) Lp = Lp + gdot_twin(i)*prm%Schmid_twin(1:3,1:3,i)
forall (k=1_pInt:3_pInt,l=1_pInt:3_pInt,m=1_pInt:3_pInt,n=1_pInt:3_pInt) & forall (k=1_pInt:3_pInt,l=1_pInt:3_pInt,m=1_pInt:3_pInt,n=1_pInt:3_pInt) &
@ -487,9 +487,9 @@ subroutine plastic_phenopowerlaw_dotState(Mp,instance,of)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! shear rates ! shear rates
call kinetics_slip(prm,stt,of,Mp,gdot_slip_pos,gdot_slip_neg) call kinetics_slip(Mp,instance,of,gdot_slip_pos,gdot_slip_neg)
dot%gamma_slip(:,of) = abs(gdot_slip_pos+gdot_slip_neg) dot%gamma_slip(:,of) = abs(gdot_slip_pos+gdot_slip_neg)
call kinetics_twin(prm,stt,of,Mp,dot%gamma_twin(:,of)) call kinetics_twin(Mp,instance,of,dot%gamma_twin(:,of))
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! hardening ! hardening
@ -509,134 +509,6 @@ subroutine plastic_phenopowerlaw_dotState(Mp,instance,of)
end subroutine plastic_phenopowerlaw_dotState end subroutine plastic_phenopowerlaw_dotState
!--------------------------------------------------------------------------------------------------
!> @brief calculates shear rates on slip systems and derivatives with respect to resolved stress
!> @details Shear rates are calculated only optionally.
! NOTE: Against the common convention, the result (i.e. intent(out)) variables are the last to
! have the optional arguments at the end
!--------------------------------------------------------------------------------------------------
pure subroutine kinetics_slip(prm,stt,of,Mp,gdot_slip_pos,gdot_slip_neg, &
dgdot_dtau_slip_pos,dgdot_dtau_slip_neg)
use prec, only: &
dNeq0
use math, only: &
math_mul33xx33
implicit none
type(tParameters), intent(in) :: &
prm
type(tPhenopowerlawState), intent(in) :: &
stt
integer(pInt), intent(in) :: &
of
real(pReal), dimension(prm%totalNslip), intent(out) :: &
gdot_slip_pos, &
gdot_slip_neg
real(pReal), dimension(prm%totalNslip), optional, intent(out) :: &
dgdot_dtau_slip_pos, &
dgdot_dtau_slip_neg
real(pReal), dimension(3,3), intent(in) :: &
Mp
real(pReal), dimension(prm%totalNslip) :: &
tau_slip_pos, &
tau_slip_neg
integer(pInt) :: i
logical :: nonSchmidActive
nonSchmidActive = size(prm%nonSchmidCoeff) > 0_pInt
do i = 1_pInt, prm%totalNslip
tau_slip_pos(i) = math_mul33xx33(Mp,prm%nonSchmid_pos(1:3,1:3,i))
tau_slip_neg(i) = merge(math_mul33xx33(Mp,prm%nonSchmid_neg(1:3,1:3,i)), &
0.0_pReal, nonSchmidActive)
enddo
where(dNeq0(tau_slip_pos))
gdot_slip_pos = prm%gdot0_slip * merge(0.5_pReal,1.0_pReal, nonSchmidActive) & ! 1/2 if non-Schmid active
* sign(abs(tau_slip_pos/stt%xi_slip(:,of))**prm%n_slip, tau_slip_pos)
else where
gdot_slip_pos = 0.0_pReal
end where
where(dNeq0(tau_slip_neg))
gdot_slip_neg = 0.5_pReal*prm%gdot0_slip &
* sign(abs(tau_slip_neg/stt%xi_slip(:,of))**prm%n_slip, tau_slip_neg)
else where
gdot_slip_neg = 0.0_pReal
end where
if (present(dgdot_dtau_slip_pos)) then
where(dNeq0(gdot_slip_pos))
dgdot_dtau_slip_pos = gdot_slip_pos*prm%n_slip/tau_slip_pos
else where
dgdot_dtau_slip_pos = 0.0_pReal
end where
endif
if (present(dgdot_dtau_slip_neg)) then
where(dNeq0(gdot_slip_neg))
dgdot_dtau_slip_neg = gdot_slip_neg*prm%n_slip/tau_slip_neg
else where
dgdot_dtau_slip_neg = 0.0_pReal
end where
endif
end subroutine kinetics_slip
!--------------------------------------------------------------------------------------------------
!> @brief calculates shear rates on twin systems and derivatives with respect to resolved stress.
! twinning is assumed to take place only in untwinned volume.
!> @details Derivates are calculated only optionally.
! NOTE: Against the common convention, the result (i.e. intent(out)) variables are the last to
! have the optional arguments at the end.
!--------------------------------------------------------------------------------------------------
pure subroutine kinetics_twin(prm,stt,of,Mp,gdot_twin,dgdot_dtau_twin)
use prec, only: &
dNeq0
use math, only: &
math_mul33xx33
implicit none
type(tParameters), intent(in) :: &
prm
type(tPhenopowerlawState), intent(in) :: &
stt
integer(pInt), intent(in) :: &
of
real(pReal), dimension(3,3), intent(in) :: &
Mp
real(pReal), dimension(prm%totalNtwin), intent(out) :: &
gdot_twin
real(pReal), dimension(prm%totalNtwin), optional, intent(out) :: &
dgdot_dtau_twin
real(pReal), dimension(prm%totalNtwin) :: &
tau_twin
integer(pInt) :: i
do i = 1_pInt, prm%totalNtwin
tau_twin(i) = math_mul33xx33(Mp,prm%Schmid_twin(1:3,1:3,i))
enddo
where(tau_twin > 0.0_pReal)
gdot_twin = (1.0_pReal-sum(stt%gamma_twin(:,of)/prm%gamma_twin_char)) & ! only twin in untwinned volume fraction
* prm%gdot0_twin*(abs(tau_twin)/stt%xi_twin(:,of))**prm%n_twin
else where
gdot_twin = 0.0_pReal
end where
if (present(dgdot_dtau_twin)) then
where(dNeq0(gdot_twin))
dgdot_dtau_twin = gdot_twin*prm%n_twin/tau_twin
else where
dgdot_dtau_twin = 0.0_pReal
end where
endif
end subroutine kinetics_twin
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief return array of constitutive results !> @brief return array of constitutive results
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
@ -674,7 +546,7 @@ function plastic_phenopowerlaw_postResults(Mp,instance,of) result(postResults)
postResults(c+1_pInt:c+prm%totalNslip) = stt%gamma_slip(1:prm%totalNslip,of) postResults(c+1_pInt:c+prm%totalNslip) = stt%gamma_slip(1:prm%totalNslip,of)
c = c + prm%totalNslip c = c + prm%totalNslip
case (shearrate_slip_ID) case (shearrate_slip_ID)
call kinetics_slip(prm,stt,of,Mp,gdot_slip_pos,gdot_slip_neg) call kinetics_slip(Mp,instance,of,gdot_slip_pos,gdot_slip_neg)
postResults(c+1_pInt:c+prm%totalNslip) = gdot_slip_pos+gdot_slip_neg postResults(c+1_pInt:c+prm%totalNslip) = gdot_slip_pos+gdot_slip_neg
c = c + prm%totalNslip c = c + prm%totalNslip
case (resolvedstress_slip_ID) case (resolvedstress_slip_ID)
@ -690,7 +562,7 @@ function plastic_phenopowerlaw_postResults(Mp,instance,of) result(postResults)
postResults(c+1_pInt:c+prm%totalNtwin) = stt%gamma_twin(1:prm%totalNtwin,of) postResults(c+1_pInt:c+prm%totalNtwin) = stt%gamma_twin(1:prm%totalNtwin,of)
c = c + prm%totalNtwin c = c + prm%totalNtwin
case (shearrate_twin_ID) case (shearrate_twin_ID)
call kinetics_twin(prm,stt,of,Mp,postResults(c+1_pInt:c+prm%totalNtwin)) call kinetics_twin(Mp,instance,of,postResults(c+1_pInt:c+prm%totalNtwin))
c = c + prm%totalNtwin c = c + prm%totalNtwin
case (resolvedstress_twin_ID) case (resolvedstress_twin_ID)
do i = 1_pInt, prm%totalNtwin do i = 1_pInt, prm%totalNtwin
@ -705,4 +577,136 @@ function plastic_phenopowerlaw_postResults(Mp,instance,of) result(postResults)
end function plastic_phenopowerlaw_postResults end function plastic_phenopowerlaw_postResults
!--------------------------------------------------------------------------------------------------
!> @brief calculates shear rates on slip systems and derivatives with respect to resolved stress
!> @details Shear rates are calculated only optionally.
! NOTE: Against the common convention, the result (i.e. intent(out)) variables are the last to
! have the optional arguments at the end
!--------------------------------------------------------------------------------------------------
pure subroutine kinetics_slip(Mp,instance,of, &
gdot_slip_pos,gdot_slip_neg,dgdot_dtau_slip_pos,dgdot_dtau_slip_neg)
use prec, only: &
dNeq0
use math, only: &
math_mul33xx33
implicit none
real(pReal), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
integer(pInt), intent(in) :: &
instance, &
of
real(pReal), dimension(param(instance)%totalNslip), intent(out) :: &
gdot_slip_pos, &
gdot_slip_neg
real(pReal), dimension(param(instance)%totalNslip), intent(out), optional :: &
dgdot_dtau_slip_pos, &
dgdot_dtau_slip_neg
real(pReal), dimension(param(instance)%totalNslip) :: &
tau_slip_pos, &
tau_slip_neg
integer(pInt) :: i
logical :: nonSchmidActive
associate(prm => param(instance), stt => state(instance))
nonSchmidActive = size(prm%nonSchmidCoeff) > 0_pInt
do i = 1_pInt, prm%totalNslip
tau_slip_pos(i) = math_mul33xx33(Mp,prm%nonSchmid_pos(1:3,1:3,i))
tau_slip_neg(i) = merge(math_mul33xx33(Mp,prm%nonSchmid_neg(1:3,1:3,i)), &
0.0_pReal, nonSchmidActive)
enddo
where(dNeq0(tau_slip_pos))
gdot_slip_pos = prm%gdot0_slip * merge(0.5_pReal,1.0_pReal, nonSchmidActive) & ! 1/2 if non-Schmid active
* sign(abs(tau_slip_pos/stt%xi_slip(:,of))**prm%n_slip, tau_slip_pos)
else where
gdot_slip_pos = 0.0_pReal
end where
where(dNeq0(tau_slip_neg))
gdot_slip_neg = 0.5_pReal*prm%gdot0_slip &
* sign(abs(tau_slip_neg/stt%xi_slip(:,of))**prm%n_slip, tau_slip_neg)
else where
gdot_slip_neg = 0.0_pReal
end where
if (present(dgdot_dtau_slip_pos)) then
where(dNeq0(gdot_slip_pos))
dgdot_dtau_slip_pos = gdot_slip_pos*prm%n_slip/tau_slip_pos
else where
dgdot_dtau_slip_pos = 0.0_pReal
end where
endif
if (present(dgdot_dtau_slip_neg)) then
where(dNeq0(gdot_slip_neg))
dgdot_dtau_slip_neg = gdot_slip_neg*prm%n_slip/tau_slip_neg
else where
dgdot_dtau_slip_neg = 0.0_pReal
end where
endif
end associate
end subroutine kinetics_slip
!--------------------------------------------------------------------------------------------------
!> @brief calculates shear rates on twin systems and derivatives with respect to resolved stress.
! twinning is assumed to take place only in untwinned volume.
!> @details Derivates are calculated only optionally.
! NOTE: Against the common convention, the result (i.e. intent(out)) variables are the last to
! have the optional arguments at the end.
!--------------------------------------------------------------------------------------------------
pure subroutine kinetics_twin(Mp,instance,of,&
gdot_twin,dgdot_dtau_twin)
use prec, only: &
dNeq0
use math, only: &
math_mul33xx33
implicit none
real(pReal), dimension(3,3), intent(in) :: &
Mp !< Mandel stress
integer(pInt), intent(in) :: &
instance, &
of
real(pReal), dimension(param(instance)%totalNtwin), intent(out) :: &
gdot_twin
real(pReal), dimension(param(instance)%totalNtwin), intent(out), optional :: &
dgdot_dtau_twin
real(pReal), dimension(param(instance)%totalNtwin) :: &
tau_twin
integer(pInt) :: i
associate(prm => param(instance), stt => state(instance))
do i = 1_pInt, prm%totalNtwin
tau_twin(i) = math_mul33xx33(Mp,prm%Schmid_twin(1:3,1:3,i))
enddo
where(tau_twin > 0.0_pReal)
gdot_twin = (1.0_pReal-sum(stt%gamma_twin(:,of)/prm%gamma_twin_char)) & ! only twin in untwinned volume fraction
* prm%gdot0_twin*(abs(tau_twin)/stt%xi_twin(:,of))**prm%n_twin
else where
gdot_twin = 0.0_pReal
end where
if (present(dgdot_dtau_twin)) then
where(dNeq0(gdot_twin))
dgdot_dtau_twin = gdot_twin*prm%n_twin/tau_twin
else where
dgdot_dtau_twin = 0.0_pReal
end where
endif
end associate
end subroutine kinetics_twin
end module plastic_phenopowerlaw end module plastic_phenopowerlaw