new formular for kappa is done, time for debugging

This commit is contained in:
zhangc43 2016-04-13 16:38:22 -04:00
parent de6b712b09
commit 00abdc34c1
3 changed files with 255 additions and 150 deletions

View File

@ -436,7 +436,7 @@ end function constitutive_homogenizedC
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief calls microstructure function of the different constitutive models !> @brief calls microstructure function of the different constitutive models
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
subroutine constitutive_microstructure(orientations, Fe, Fp, ipc, ip, el, F0s,Fes,Fps) subroutine constitutive_microstructure(orientations, Fe, Fp, ipc, ip, el, F0s,Fes,Fps,Tstar_vs)
use prec, only: & use prec, only: &
pReal pReal
use material, only: & use material, only: &
@ -473,10 +473,15 @@ subroutine constitutive_microstructure(orientations, Fe, Fp, ipc, ip, el, F0s,Fe
ho, & !< homogenization ho, & !< homogenization
tme !< thermal member position tme !< thermal member position
real(pReal), intent(in), dimension(:,:,:,:) :: & real(pReal), intent(in), dimension(:,:,:,:) :: &
orientations, & orientations
real(pReal), intent(in), dimension(:,:,:,:,:) :: &
F0s, & F0s, &
Fes, & Fes, &
Fps !< crystal orientations as quaternions Fps
real(pReal), intent(in), dimension(6,:,:,:) :: &
Tstar_vs !< crystal orientations as quaternions
ho = material_homog(ip,el) ho = material_homog(ip,el)
tme = thermalMapping(ho)%p(ip,el) tme = thermalMapping(ho)%p(ip,el)
@ -491,7 +496,7 @@ subroutine constitutive_microstructure(orientations, Fe, Fp, ipc, ip, el, F0s,Fe
case (PLASTICITY_NONLOCAL_ID) plasticityType case (PLASTICITY_NONLOCAL_ID) plasticityType
call plastic_nonlocal_microstructure (Fe,Fp,ip,el) call plastic_nonlocal_microstructure (Fe,Fp,ip,el)
case (PLASTICITY_PHENOPLUS_ID) plasticityType case (PLASTICITY_PHENOPLUS_ID) plasticityType
call plastic_phenoplus_microstructure(orientations,ipc,ip,el,F0s,Fes,Fps) call plastic_phenoplus_microstructure(orientations,ipc,ip,el,F0s,Fes,Fps,Tstar_vs)
end select plasticityType end select plasticityType
end subroutine constitutive_microstructure end subroutine constitutive_microstructure

View File

@ -440,7 +440,8 @@ subroutine crystallite_init
c,i,e, c,i,e,
crystallite_F0, crystallite_F0,
crystallite_Fe, crystallite_Fe,
crystallite_Fp) ! update dependent state variables to be consistent with basic states crystallite_Fp,
crystallite_Tstar_v) ! update dependent state variables to be consistent with basic states
enddo enddo
enddo enddo
enddo enddo
@ -1720,7 +1721,8 @@ subroutine crystallite_integrateStateRK4()
g, i, e, g, i, e,
crystallite_F0, crystallite_F0,
crystallite_Fe, crystallite_Fe,
crystallite_Fp) ! update dependent state variables to be consistent with basic states crystallite_Fp,
crystallite_Tstar_v) ! update dependent state variables to be consistent with basic states
enddo; enddo; enddo enddo; enddo; enddo
!$OMP ENDDO !$OMP ENDDO
@ -2049,7 +2051,8 @@ subroutine crystallite_integrateStateRKCK45()
g, i, e, g, i, e,
crystallite_F0, crystallite_F0,
crystallite_Fe, crystallite_Fe,
crystallite_Fp) ! update dependent state variables to be consistent with basic states crystallite_Fp,
crystallite_Tstar_v) ! update dependent state variables to be consistent with basic states
enddo; enddo; enddo enddo; enddo; enddo
!$OMP ENDDO !$OMP ENDDO
@ -2272,7 +2275,8 @@ subroutine crystallite_integrateStateRKCK45()
g, i, e, g, i, e,
crystallite_F0, crystallite_F0,
crystallite_Fe, crystallite_Fe,
crystallite_Fp) ! update dependent state variables to be consistent with basic states crystallite_Fp,
crystallite_Tstar_v) ! update dependent state variables to be consistent with basic states
enddo; enddo; enddo enddo; enddo; enddo
!$OMP ENDDO !$OMP ENDDO
@ -2510,7 +2514,8 @@ subroutine crystallite_integrateStateAdaptiveEuler()
g, i, e, g, i, e,
crystallite_F0, crystallite_F0,
crystallite_Fe, crystallite_Fe,
crystallite_Fp) ! update dependent state variables to be consistent with basic states crystallite_Fp,
crystallite_Tstar_v) ! update dependent state variables to be consistent with basic states
enddo; enddo; enddo enddo; enddo; enddo
!$OMP ENDDO !$OMP ENDDO
!$OMP END PARALLEL !$OMP END PARALLEL
@ -2857,7 +2862,8 @@ eIter = FEsolving_execElem(1:2)
g, i, e, g, i, e,
crystallite_F0, crystallite_F0,
crystallite_Fe, crystallite_Fe,
crystallite_Fp) ! update dependent state variables to be consistent with basic states crystallite_Fp,
crystallite_Tstar_v) ! update dependent state variables to be consistent with basic states
enddo; enddo; enddo enddo; enddo; enddo
!$OMP ENDDO !$OMP ENDDO
!$OMP END PARALLEL !$OMP END PARALLEL
@ -3105,7 +3111,8 @@ subroutine crystallite_integrateStateFPI()
g, i, e, g, i, e,
crystallite_F0, crystallite_F0,
crystallite_Fe, crystallite_Fe,
crystallite_Fp) ! update dependent state variables to be consistent with basic states crystallite_Fp,
crystallite_Tstar_v) ! update dependent state variables to be consistent with basic states
p = phaseAt(g,i,e) p = phaseAt(g,i,e)
c = phasememberAt(g,i,e) c = phasememberAt(g,i,e)
plasticState(p)%previousDotState2(:,c) = plasticState(p)%previousDotState(:,c) plasticState(p)%previousDotState2(:,c) = plasticState(p)%previousDotState(:,c)

View File

@ -739,162 +739,255 @@ end subroutine plastic_phenoplus_aTolState
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief calculate push-up factors (kappa) for each voxel based on its neighbors !> @brief calculate push-up factors (kappa) for each voxel based on its neighbors
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
subroutine plastic_phenoplus_microstructure(orientation,ipc,ip,el,F0,Fe,Fp) subroutine plastic_phenoplus_microstructure(orientation,ipc,ip,el,F0,Fe,Fp,Tstar_v)
use math, only: pi, & use math, only: pi, &
math_identity2nd, & math_identity2nd, &
math_mul33x33, & math_mul33x33, &
math_mul33xx33, & math_mul33xx33, &
math_mul3x3, & math_mul3x3, &
math_transpose33, & math_transpose33, &
math_qDot, & math_qDot, &
math_qRot, & math_qRot, &
indeg indeg
use mesh, only: mesh_element, & use mesh, only: mesh_element, &
FE_NipNeighbors, & FE_NipNeighbors, &
FE_geomtype, & FE_geomtype, &
FE_celltype, & FE_celltype, &
mesh_maxNips, & mesh_maxNips, &
mesh_NcpElems, & mesh_NcpElems, &
mesh_ipNeighborhood mesh_ipNeighborhood
use material, only: material_phase, & use material, only: material_phase, &
material_texture, & material_texture, &
phase_plasticityInstance, & phase_plasticityInstance, &
phaseAt, phasememberAt, & phaseAt, phasememberAt, &
homogenization_maxNgrains, & homogenization_maxNgrains, &
plasticState plasticState
use lattice, only: lattice_sn, & use lattice, only: lattice_Sslip_v, &
lattice_sd, & lattice_maxNslipFamily, &
lattice_qDisorientation lattice_NslipSystem, &
lattice_NslipSystem, &
lattice_sn, &
lattice_sd, &
lattice_qDisorientation
!***input variables !***input variables
implicit none implicit none
integer(pInt), intent(in) :: & integer(pInt), intent(in) :: &
ipc, & !< component-ID of integration point ipc, & !< component-ID of integration point
ip, & !< integration point ip, & !< integration point
el el
real(pReal), dimension(3,3,homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems), intent(in) :: & real(pReal), dimension(3,3,homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems), intent(in) :: &
F0, & ! deformation gradient from last increment F0, & !< deformation gradient from last increment
Fe, & ! elastic deformation gradient Fe, & !< elastic deformation gradient
Fp ! elastic deformation gradient !< element Fp !< elastic deformation gradient !< element
real(pReal), dimension(4,homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems), intent(in) :: & real(pReal), dimension(4,homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems), intent(in) :: &
orientation ! crystal orientation in quaternions orientation !< crystal orientation in quaternions
real(pReal), dimension(6,homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems), intent(in) :: &
Tstar_v !< for calculation of gdot
!***local variables !***local variables
integer(pInt) instance, & !my instance of this plasticity integer(pInt) instance, & !my instance of this plasticity
ph, & !my phase ph, & !my phase
of, & !my spatial position in memory (offset) of, & !my spatial position in memory (offset)
textureID, & !my texture textureID, & !my texture
Nneighbors, & !number of neighbors (<= 6) Nneighbors, & !number of neighbors (<= 6)
vld_Nneighbors, & !number of my valid neighbors vld_Nneighbors, & !number of my valid neighbors
n, & !neighbor index (for iterating through all neighbors) n, & !neighbor index (for iterating through all neighbors)
ns, & !number of slip system n_calcTaylor, & !
nt, & !number of twin system n_phasecheck, & !
me_slip, & !my slip system index ns, & !number of slip system
neighbor_el, & !element number of neighboring material point nt, & !number of twin system
neighbor_ip, & !integration point of neighboring material point me_slip, & !my slip system index
neighbor_n, & !I have no idea what is this neighbor_el, & !element number of neighboring material point
neighbor_of, & !spatial position in memory for this neighbor (offset) neighbor_ip, & !integration point of neighboring material point
neighbor_ph, & !neighbor's phase neighbor_ipc, & !I have no idea what is this
neighbor_tex, & !neighbor's texture ID neighbor_of, & !spatial position in memory for this neighbor (offset)
ne_slip_ac, & !loop to find neighbor shear neighbor_ph, & !neighbor's phase
ne_slip, & !slip system index for neighbor neighbor_instance, & !neighbor's instance of this plasticity
index_kappa, & !index of pushup factors in plasticState neighbor_tex, & !neighbor's texture ID
offset_acshear_slip, & !offset in PlasticState for the accumulative shear ne_slip, & !slip system index for neighbor
j !quickly loop through slip families index_kappa, & !index of pushup factors in plasticState
j, & !quickly loop through slip families
f,i,& !loop counter for me
f_ne, i_ne !loop counter for neighbor
real(pReal) kappa_max, & ! real(pReal) kappa_max, & !
tmp_myshear_slip, & !temp storage for accumulative shear for me tmp_myshear_slip, & !temp storage for accumulative shear for me
mprime_cut, & !m' cutoff to consider neighboring effect mprime_cut, & !m' cutoff to consider neighboring effect
dtaylor_cut, & !threshold for determine high contrast interface using Taylor factor dtaylor_cut, & !threshold for determine high contrast interface using Taylor factor
avg_acshear_ne, & !the average accumulative shear from my neighbor avg_acshear_ne, & !the average accumulative shear from my neighbor
taylor_me, & !Taylor factor for me taylor_me, & !Taylor factor for me
taylor_ne, & !Taylor factor for my current neighbor taylor_ne, & !Taylor factor for my current neighbor
tmp_mprime, & !temp holder for m' value d_vonstrain !von Mises delta strain (temp container)
tmp_acshear !temp holder for accumulative shear for m'
real(pReal), dimension(3,3) :: & real(pReal), dimension(3,3) :: &
F0_me, & !my deformation gradient from last converged increment F0_me, & !my deformation gradient from last converged increment
Fe_me, & !my elastic deformation gradient Fe_me, & !my elastic deformation gradient
Fp_me, & !my plastic deformation gradient Fp_me, & !my plastic deformation gradient
dF_me, & !my deformation gradient change (delta) dF_me, & !my deformation gradient change (delta)
dE_me, & !my Green Lagrangian strain tensor (delta) dE_me, & !my Green Lagrangian strain tensor (delta)
Fe_ne, & !elastic deformation gradient of my current neighbor Fe_ne, & !elastic deformation gradient of my current neighbor
Fp_ne, & !plastic deformation gradient of my current neighbor Fp_ne, & !plastic deformation gradient of my current neighbor
dF_ne, & !deformation gradient of my current neighbor dF_ne, & !deformation gradient of my current neighbor
dE_ne !delta Green Lagrangian strain tensor dE_ne !delta Green Lagrangian strain tensor
real(pReal), dimension(plastic_phenoplus_totalNslip(phase_plasticityInstance(material_phase(1,ip,el)))) :: & real(pReal), dimension(plastic_phenoplus_totalNslip(phase_plasticityInstance(material_phase(1,ip,el)))) :: &
m_primes, & !m' between me_alpha(one) and neighbor beta(all) m_primes !m' between me_alpha(one) and neighbor beta(all)
me_acshear, & !temp storage for ac_shear of one particular system for me
ne_acshear !temp storage for ac_shear of one particular system for one of my neighbor
real(pReal), dimension(3,plastic_phenoplus_totalNslip(phase_plasticityInstance(material_phase(1,ip,el)))) :: & real(pReal), dimension(3,plastic_phenoplus_totalNslip(phase_plasticityInstance(material_phase(1,ip,el)))) :: &
slipNormal, & slipNormal, &
slipDirect slipDirect
real(pReal), dimension(4) :: my_orientation, & !store my orientation real(pReal), dimension(4) :: &
neighbor_orientation, & !store my neighbor orientation my_orientation, & !store my orientation
absMisorientation neighbor_orientation, & !store my neighbor orientation
absMisorientation
real(pReal), dimension(FE_NipNeighbors(FE_celltype(FE_geomtype(mesh_element(2,el))))) :: & real(pReal), dimension(FE_NipNeighbors(FE_celltype(FE_geomtype(mesh_element(2,el))))) :: &
ne_mprimes !m' between each neighbor ne_mprimes, & !m' between each neighbor
d_taylors !store (taylor_ne-taylor_me) for each neighbor
!***Get my properties !***Get my properties
Nneighbors = FE_NipNeighbors(FE_celltype(FE_geomtype(mesh_element(2,el)))) Nneighbors = FE_NipNeighbors(FE_celltype(FE_geomtype(mesh_element(2,el))))
ph = phaseAt(ipc,ip,el) !get my phase ph = phaseAt(ipc,ip,el) !get my phase
of = phasememberAt(ipc,ip,el) !get my spatial location offset in memory of = phasememberAt(ipc,ip,el) !get my spatial location offset in memory
textureID = material_texture(1,ip,el) !get my texture ID textureID = material_texture(1,ip,el) !get my texture ID
instance = phase_plasticityInstance(ph) !get my instance based on phase ID instance = phase_plasticityInstance(ph) !get my instance based on phase ID
ns = plastic_phenoplus_totalNslip(instance) ns = plastic_phenoplus_totalNslip(instance)
nt = plastic_phenoplus_totalNtwin(instance) nt = plastic_phenoplus_totalNtwin(instance)
offset_acshear_slip = ns + nt + 2_pInt index_kappa = ns + nt + 2_pInt + ns + nt !location of kappa in plasticState
index_kappa = ns + nt + 2_pInt + ns + nt !location of kappa in plasticState
!***init calculation for given voxel !***init calculation for given voxel
mprime_cut = 0.7_pReal !set by Dr.Bieler mprime_cut = 0.7_pReal !set by Dr.Bieler
dtaylor_cut = 1.0_pReal !set by Chen, quick test only dtaylor_cut = 1.0_pReal !set by Chen, quick test only
!***gather my orientation, F and slip systems !***gather my orientation, F and slip systems
my_orientation = orientation(1:4, ipc, ip, el) my_orientation = orientation(1:4, ipc, ip, el)
F0_me = F0(1:3, 1:3, ipc, ip, el) F0_me = F0(1:3, 1:3, ipc, ip, el)
Fe_me = Fe(1:3, 1:3, ipc, ip, el) Fe_me = Fe(1:3, 1:3, ipc, ip, el)
Fp_me = Fp(1:3, 1:3, ipc, ip, el) Fp_me = Fp(1:3, 1:3, ipc, ip, el)
slipNormal(1:3, 1:ns) = lattice_sn(1:3, 1:ns, ph) slipNormal(1:3, 1:ns) = lattice_sn(1:3, 1:ns, ph)
slipDirect(1:3, 1:ns) = lattice_sd(1:3, 1:ns, ph) slipDirect(1:3, 1:ns) = lattice_sd(1:3, 1:ns, ph)
!******calculate Taylor factor for me
!@note: we need teh !***check if all my neighbors have the same phase as me
F_me = math_mul33x33(Fe_me,Fp_me) vld_Nneighbors = 0
E_me = 0.5*(math_mul33x33(math_transpose33(F_me), F_me) - math_identity2nd) !E = 0.5(F^tF-I) PHASECHECK DO n_phasecheck = 1_pInt, Nneighbors
vonStrain !******for each of my neighbor
neighbor_el = mesh_ipNeighborhood( 1, n_phasecheck, ip, el )
neighbor_ip = mesh_ipNeighborhood( 2, n_phasecheck, ip, el )
neighbor_ipc = 1
neighbor_of = phasememberAt( neighbor_ipc, neighbor_ip, neighbor_el )
neighbor_ph = phaseAt( neighbor_ipc, neighbor_ip, neighbor_el )
IF (neighbor_ph == ph) THEN
vld_Nneighbors = vld_Nneighbors + 1_pInt
ENDIF
ENDDO PHASECHECK
!***initialize kappa with 1.0 (assume no push-up)
plasticState(ph)%state(index_kappa+1_pInt:index_kappa+ns, of) = 1.0_pReal
!***only calculate kappa for those inside the main phase
IF (vld_Nneighbors == Nneighbors) THEN
!******calculate Taylor factor for me
dF_me = math_mul33x33(Fe_me,Fp_me) - F0_me
dE_me = 0.5*(math_mul33x33(math_transpose33(dF_me), dF_me) - math_identity2nd(3)) !dE = 0.5(dF^tdF-I)
d_vonstrain = SQRT(2.0_pReal/3.0_pReal * math_mul33xx33(dE_me, dE_me))
sum_gdot = 0.0_pReal
!go through my slip system to find the sum of gamma_dot
j = 0_pInt
slipFamilies: DO f = 1_pInt,lattice_maxNslipFamily
index_myFamily = sum(lattice_NslipSystem(1:f-1_pInt,ph)) !at which index starts my family
slipSystems: DO i = 1_pInt,plastic_phenoplus_Nslip(f,instance)
j = j+1_pInt
tau_slip = dot_product(Tstar_v(1:6, ipc, ip, el),lattice_Sslip_v(1:6,1,index_myFamily+i,ph))
sum_gdot = sum_gdot + &
plastic_phenoplus_gdot0_slip(instance)* &
((abs(tau_slip)/(state(instance)%s_slip(j,of))) &
**plastic_phenoplus_n_slip(instance))*sign(1.0_pReal,tau_slip)
ENDDO slipSystems
ENDDO slipFamilies
taylor_me = d_vonstrain/sum_gdot
!***calculate delta_M (Taylor factor) between each neighbor and me
LOOPCALCTAYLOR: DO n_calcTaylor=1_pInt, Nneighbors
!******for each of my neighbor
neighbor_el = mesh_ipNeighborhood( 1, n_calcTaylor, ip, el )
neighbor_ip = mesh_ipNeighborhood( 2, n_calcTaylor, ip, el )
neighbor_ipc = 1 !It is ipc
neighbor_of = phasememberAt( neighbor_ipc, neighbor_ip, neighbor_el )
neighbor_ph = phaseAt( neighbor_ipc, neighbor_ip, neighbor_el )
neighbor_instance = phase_plasticityInstance( neighbor_ph )
neighbor_tex = material_texture( 1,neighbor_ip, neighbor_el )
neighbor_orientation = orientation( 1:4, neighbor_ipc, neighbor_ip, neighbor_el ) !ipc is always 1.
Fe_ne = Fe( 1:3, 1:3, neighbor_ipc, neighbor_ip, neighbor_el )
Fp_ne = Fp( 1:3, 1:3, neighbor_ipc, neighbor_ip, neighbor_el )
F0_ne = F0( 1:3, 1:3, neighbor_ipc, neighbor_ip, neighbor_el )
!******calculate the Taylor factor
dF_ne = math_mul33x33(Fe_ne, Fp_ne) - F0_ne
dE_ne = 0.5*(math_mul33x33(math_transpose33(dF_ne), dF_ne) - math_identity2nd(3)) !dE = 0.5(dF^tdF-I)
d_vonstrain = SQRT(2.0_pReal/3.0_pReal * math_mul33xx33(dE_ne, dE_ne))
sum_gdot = 0.0_pReal
!go through my neighbor slip system to calculate sum_gdot
j = 0_pInt
slipFamiliesNeighbor: DO f_ne = 1_pInt,lattice_maxNslipFamily
index_myFamily = sum(lattice_NslipSystem(1:f_ne-1_pInt,neighbor_ph)) ! at which index starts my family
slipSystemsNeighbor: DO i_ne = 1_pInt,plastic_phenopowerlaw_Nslip(f_ne,neighbor_instance)
j = j+1_pInt
tau_slip = dot_product(Tstar_v(1:6, neighbor_ipc, neighbor_ip, neighbor_el),
lattice_Sslip_v(1:6,1,index_myFamily+i_ne,neighbor_ph))
sum_gdot = sum_gdot &
+plastic_phenopowerlaw_gdot0_slip(neighbor_instance) &
*((abs(tau_slip)/(state(neighbor_instance)%s_slip(j,neighbor_of))) &
**plastic_phenopowerlaw_n_slip(neighbor_instance))*sign(1.0_pReal,tau_slip)
ENDDO slipSystemsNeighbor
ENDDO slipFamiliesNeighbor
taylor_ne = d_vonstrain / sum_gdot
!******calculate Taylor difference
d_taylors(n_calcTaylor) = taylor_ne - taylor_me
ENDDO LOOPCALCTAYLOR
!***Only perform necessary calculation if high contrast interface is detected
IF (max(d_taylors) > dtaylor_cut) THEN
!*****calculate kappa per slip system base
LOOPMYSLIP DO me_slip = 1_pInt, ns
ne_mprimes = 0.0_pReal !initialize max m' to 0 for all neighbors
LOOPMYNEIGHBORS DO n=1_pInt, Nneighbors
!*******only consider neighbor at the high contrast interface
IF (d_taylors(n) > dtaylor_cut) THEN
neighbor_el = mesh_ipNeighborhood( 1, n_calcTaylor, ip, el )
neighbor_ip = mesh_ipNeighborhood( 2, n_calcTaylor, ip, el )
neighbor_ipc = 1 !It is ipc
neighbor_of = phasememberAt( neighbor_ipc, neighbor_ip, neighbor_el )
neighbor_ph = phaseAt( neighbor_ipc, neighbor_ip, neighbor_el )
neighbor_instance = phase_plasticityInstance( neighbor_ph )
neighbor_tex = material_texture( 1,neighbor_ip, neighbor_el )
neighbor_orientation = orientation( 1:4, neighbor_ipc, neighbor_ip, neighbor_el ) !ipc is always 1.
absMisorientation = lattice_qDisorientation( my_orientation, &
neighbor_orientation, &
0_pInt ) !no need for explicit calculation of symmetry
!*********go through neighbor slip system to calculate m'
LOOPNEIGHBORSLIP: DO ne_slip=1_pInt,ns
m_primes(ne_slip) = abs(math_mul3x3(slipNormal(1:3,me_slip), &
math_qRot(absMisorientation, slipNormal(1:3,ne_slip)))) &
*abs(math_mul3x3(slipDirect(1:3,me_slip), &
math_qRot(absMisorientation, slipDirect(1:3,ne_slip))))
ENDDO LOOPNEIGHBORSLIP
ne_mprimes(n) = max(m_primes)
ENDIF
!*******check if one of the neighbor already can provide a kick for this slip system
IF ( max(ne_mprimes) > mprime_cut ) THEN
plasticState(ph)%state(index_kappa+me_slip, of) = 1.5_pReal
EXIT
ENDIF
ENDDO LOOPMYNEIGHBORS
ENDDO LOOPMYSLIP
ENDIF
!***loop into the geometry to figure out who is my closest neighbor
LOOPNEIGHBORS: DO n=1_pInt, Nneighbors
!******for each of my neighbor, calculate the Taylor factor
ne_taylor = 1.0
!*********for the high contrast interface
IF (abs(taylor_ne - taylor_me) > dtaylor_cut) THEN
!********* gather neighbor orientation and slip systems
!********* calculate m' (need to loop through all my slip systems as well)
!********* if m'>mprime_cut kappa=1.5 else 1.0
!******
ELSE
ENDIF ENDIF
!***end of search
ENDDO LOOPNEIGHBORS
! !***gather my accumulative shear from palsticState
! FINDMYSHEAR: do j = 1_pInt,ns
! me_acshear(j) = plasticState(ph)%state(offset_acshear_slip+j, of)
! enddo FINDMYSHEAR
! !***gather my orientation and slip systems
! my_orientation = orientation(1:4, ipc, ip, el)
! slipNormal(1:3, 1:ns) = lattice_sn(1:3, 1:ns, ph)
! slipDirect(1:3, 1:ns) = lattice_sd(1:3, 1:ns, ph)
! kappa_max = plastic_phenoplus_kappa_max(instance) !maximum pushups allowed (READIN)
! !***calculate kappa between me and all my neighbors ! !***calculate kappa between me and all my neighbors
! LOOPMYSLIP: DO me_slip=1_pInt,ns ! LOOPMYSLIP: DO me_slip=1_pInt,ns