1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/constitutive_plastic_kinehardening.f90"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/constitutive_plastic_kinehardening.f90"
18 real(pReal),
allocatable,
dimension(:) :: &
19 theta0, & !< initial hardening rate of forward stress for each slip
20 theta1, & !< asymptotic hardening rate of forward stress for each slip
21 theta0_b, & !< initial hardening rate of back stress for each slip
22 theta1_b, & !< asymptotic hardening rate of back stress for each slip
25 real(pReal),
allocatable,
dimension(:,:) :: &
27 real(pReal),
allocatable,
dimension(:,:,:) :: &
32 sum_N_sl, & !< total number of active slip system
35 nonSchmidActive = .false.
36 character(len=pStringLen),
allocatable,
dimension(:) :: &
40 type :: tkinehardeningstate
41 real(pReal),
pointer,
dimension(:,:) :: & !< vectors along NipcMyInstance
42 crss, & !< critical resolved stress
43 crss_back, & !< critical resolved back stress
44 sense, & !< sense of acting shear stress (-1 or +1)
45 chi0, & !< backstress at last switch of stress sense
46 gamma0, & !< accumulated shear at last switch of stress sense
48 end type tkinehardeningstate
52 type(tParameters),
allocatable,
dimension(:) :: param
53 type(tKinehardeningState),
allocatable,
dimension(:) :: &
65 module subroutine plastic_kinehardening_init
71 sizeState, sizeDeltaState, sizeDotState, &
73 integer,
dimension(:),
allocatable :: &
75 real(pReal),
dimension(:),
allocatable :: &
76 xi_0, & !< initial resistance against plastic flow
78 character(len=pStringLen) :: &
81 write(6,
'(/,a)')
' <<<+- plastic_'//plasticity_kinehardening_label//
' init -+>>>';
flush(6)
83 ninstance = count(phase_plasticity == plasticity_kinehardening_id)
84 if (iand(debug_level(debug_constitutive),debug_levelbasic) /= 0) &
85 write(6,
'(a16,1x,i5,/)')
'# instances:',ninstance
87 allocate(param(ninstance))
88 allocate(state(ninstance))
89 allocate(dotstate(ninstance))
90 allocate(deltastate(ninstance))
92 do p = 1,
size(phase_plasticityinstance)
93 if (phase_plasticity(p) /= plasticity_kinehardening_id) cycle
94 associate(prm => param(phase_plasticityinstance(p)), &
95 dot => dotstate(phase_plasticityinstance(p)), &
96 dlt => deltastate(phase_plasticityinstance(p)), &
97 stt => state(phase_plasticityinstance(p)),&
100 prm%output =
config%getStrings(
'(output)',defaultval=emptystringarray)
110 n_sl =
config%getInts(
'nslip',defaultval=emptyintarray)
111 prm%sum_N_sl = sum(abs(n_sl))
112 slipactive:
if (prm%sum_N_sl > 0)
then
113 prm%P = lattice_schmidmatrix_slip(n_sl,
config%getString(
'lattice_structure'),&
114 config%getFloat(
'c/a',defaultval=0.0_preal))
116 if(trim(
config%getString(
'lattice_structure')) ==
'bcc')
then
117 a =
config%getFloats(
'nonschmid_coefficients',defaultval = emptyrealarray)
118 if(
size(a) > 0) prm%nonSchmidActive = .true.
119 prm%nonSchmid_pos = lattice_nonschmidmatrix(n_sl,a,+1)
120 prm%nonSchmid_neg = lattice_nonschmidmatrix(n_sl,a,-1)
122 prm%nonSchmid_pos = prm%P
123 prm%nonSchmid_neg = prm%P
125 prm%interaction_SlipSlip = lattice_interaction_slipbyslip(n_sl, &
126 config%getFloats(
'interaction_slipslip'), &
127 config%getString(
'lattice_structure'))
129 xi_0 =
config%getFloats(
'crss0', requiredsize=
size(n_sl))
130 prm%tau1 =
config%getFloats(
'tau1', requiredsize=
size(n_sl))
131 prm%tau1_b =
config%getFloats(
'tau1_b', requiredsize=
size(n_sl))
132 prm%theta0 =
config%getFloats(
'theta0', requiredsize=
size(n_sl))
133 prm%theta1 =
config%getFloats(
'theta1', requiredsize=
size(n_sl))
134 prm%theta0_b =
config%getFloats(
'theta0_b', requiredsize=
size(n_sl))
135 prm%theta1_b =
config%getFloats(
'theta1_b', requiredsize=
size(n_sl))
137 prm%gdot0 =
config%getFloat(
'gdot0')
138 prm%n =
config%getFloat(
'n_slip')
141 xi_0 = math_expand(xi_0, n_sl)
142 prm%tau1 = math_expand(prm%tau1, n_sl)
143 prm%tau1_b = math_expand(prm%tau1_b, n_sl)
144 prm%theta0 = math_expand(prm%theta0, n_sl)
145 prm%theta1 = math_expand(prm%theta1, n_sl)
146 prm%theta0_b = math_expand(prm%theta0_b,n_sl)
147 prm%theta1_b = math_expand(prm%theta1_b,n_sl)
151 if ( prm%gdot0 <= 0.0_preal) extmsg = trim(extmsg)//
' gdot0'
152 if ( prm%n <= 0.0_preal) extmsg = trim(extmsg)//
' n_slip'
153 if (any(xi_0 <= 0.0_preal)) extmsg = trim(extmsg)//
' crss0'
154 if (any(prm%tau1 <= 0.0_preal)) extmsg = trim(extmsg)//
' tau1'
155 if (any(prm%tau1_b <= 0.0_preal)) extmsg = trim(extmsg)//
' tau1_b'
159 xi_0 = emptyrealarray
160 allocate(prm%tau1,prm%tau1_b,prm%theta0,prm%theta1,prm%theta0_b,prm%theta1_b,source=emptyrealarray)
161 allocate(prm%interaction_SlipSlip(0,0))
166 nipcmyphase = count(material_phaseat == p) * discretization_nip
167 sizedotstate =
size([
'crss ',
'crss_back',
'accshear ']) * prm%sum_N_sl
168 sizedeltastate =
size([
'sense ',
'chi0 ',
'gamma0' ]) * prm%sum_N_sl
169 sizestate = sizedotstate + sizedeltastate
171 call material_allocateplasticstate(p,nipcmyphase,sizestate,sizedotstate,sizedeltastate)
176 endindex = prm%sum_N_sl
177 stt%crss => plasticstate(p)%state (startindex:endindex,:)
178 stt%crss = spread(xi_0, 2, nipcmyphase)
179 dot%crss => plasticstate(p)%dotState(startindex:endindex,:)
180 plasticstate(p)%atol(startindex:endindex) =
config%getFloat(
'atol_xi',defaultval=1.0_preal)
181 if(any(plasticstate(p)%atol(startindex:endindex) < 0.0_preal)) extmsg = trim(extmsg)//
' atol_xi'
183 startindex = endindex + 1
184 endindex = endindex + prm%sum_N_sl
185 stt%crss_back => plasticstate(p)%state (startindex:endindex,:)
186 dot%crss_back => plasticstate(p)%dotState(startindex:endindex,:)
187 plasticstate(p)%atol(startindex:endindex) =
config%getFloat(
'atol_xi',defaultval=1.0_preal)
189 startindex = endindex + 1
190 endindex = endindex + prm%sum_N_sl
191 stt%accshear => plasticstate(p)%state (startindex:endindex,:)
192 dot%accshear => plasticstate(p)%dotState(startindex:endindex,:)
193 plasticstate(p)%atol(startindex:endindex) =
config%getFloat(
'atol_gamma',defaultval=1.0e-6_preal)
194 if(any(plasticstate(p)%atol(startindex:endindex) < 0.0_preal)) extmsg = trim(extmsg)//
' atol_gamma'
196 plasticstate(p)%slipRate => plasticstate(p)%dotState(startindex:endindex,:)
198 o = plasticstate(p)%offsetDeltaState
199 startindex = endindex + 1
200 endindex = endindex + prm%sum_N_sl
201 stt%sense => plasticstate(p)%state (startindex :endindex ,:)
202 dlt%sense => plasticstate(p)%deltaState(startindex-o:endindex-o,:)
204 startindex = endindex + 1
205 endindex = endindex + prm%sum_N_sl
206 stt%chi0 => plasticstate(p)%state (startindex :endindex ,:)
207 dlt%chi0 => plasticstate(p)%deltaState(startindex-o:endindex-o,:)
209 startindex = endindex + 1
210 endindex = endindex + prm%sum_N_sl
211 stt%gamma0 => plasticstate(p)%state (startindex :endindex ,:)
212 dlt%gamma0 => plasticstate(p)%deltaState(startindex-o:endindex-o,:)
214 plasticstate(p)%state0 = plasticstate(p)%state
220 if (extmsg /=
'')
call io_error(211,ext_msg=trim(extmsg)//
'('//plasticity_kinehardening_label//
')')
225 end subroutine plastic_kinehardening_init
231 pure module subroutine plastic_kinehardening_lpanditstangent(lp,dlp_dmp,mp,instance,of)
233 real(pReal),
dimension(3,3),
intent(out) :: &
235 real(pReal),
dimension(3,3,3,3),
intent(out) :: &
238 real(pReal),
dimension(3,3),
intent(in) :: &
240 integer,
intent(in) :: &
246 real(pReal),
dimension(param(instance)%sum_N_sl) :: &
248 dgdot_dtau_pos,dgdot_dtau_neg
253 associate(prm => param(instance))
255 call kinetics(mp,instance,of,gdot_pos,gdot_neg,dgdot_dtau_pos,dgdot_dtau_neg)
257 do i = 1, prm%sum_N_sl
258 lp = lp + (gdot_pos(i)+gdot_neg(i))*prm%P(1:3,1:3,i)
259 forall (k=1:3,l=1:3,m=1:3,n=1:3) &
260 dlp_dmp(k,l,m,n) = dlp_dmp(k,l,m,n) &
261 + dgdot_dtau_pos(i) * prm%P(k,l,i) * prm%nonSchmid_pos(m,n,i) &
262 + dgdot_dtau_neg(i) * prm%P(k,l,i) * prm%nonSchmid_neg(m,n,i)
267 end subroutine plastic_kinehardening_lpanditstangent
273 module subroutine plastic_kinehardening_dotstate(mp,instance,of)
275 real(pReal),
dimension(3,3),
intent(in) :: &
277 integer,
intent(in) :: &
283 real(pReal),
dimension(param(instance)%sum_N_sl) :: &
287 associate(prm => param(instance), stt => state(instance), dot => dotstate(instance))
289 call kinetics(mp,instance,of,gdot_pos,gdot_neg)
290 dot%accshear(:,of) = abs(gdot_pos+gdot_neg)
291 sumgamma = sum(stt%accshear(:,of))
294 dot%crss(:,of) = matmul(prm%interaction_SlipSlip,dot%accshear(:,of)) &
296 + (prm%theta0 - prm%theta1 + prm%theta0*prm%theta1*sumgamma/prm%tau1) &
297 * exp(-sumgamma*prm%theta0/prm%tau1) &
300 dot%crss_back(:,of) = stt%sense(:,of)*dot%accshear(:,of) * &
302 (prm%theta0_b - prm%theta1_b &
303 + prm%theta0_b*prm%theta1_b/(prm%tau1_b+stt%chi0(:,of))*(stt%accshear(:,of)-stt%gamma0(:,of))&
304 ) *exp(-(stt%accshear(:,of)-stt%gamma0(:,of)) *prm%theta0_b/(prm%tau1_b+stt%chi0(:,of))) &
309 end subroutine plastic_kinehardening_dotstate
315 module subroutine plastic_kinehardening_deltastate(mp,instance,of)
317 real(pReal),
dimension(3,3),
intent(in) :: &
319 integer,
intent(in) :: &
323 real(pReal),
dimension(param(instance)%sum_N_sl) :: &
327 associate(prm => param(instance), stt => state(instance), dlt => deltastate(instance))
329 call kinetics(mp,instance,of,gdot_pos,gdot_neg)
330 sense = merge(state(instance)%sense(:,of), &
331 sign(1.0_preal,gdot_pos+gdot_neg), &
332 deq0(gdot_pos+gdot_neg,1e-10_preal))
334 # 338 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/constitutive_plastic_kinehardening.f90"
338 where(dneq(sense,stt%sense(:,of),0.1_preal))
339 dlt%sense (:,of) = sense - stt%sense(:,of)
340 dlt%chi0 (:,of) = abs(stt%crss_back(:,of)) - stt%chi0(:,of)
341 dlt%gamma0(:,of) = stt%accshear(:,of) - stt%gamma0(:,of)
343 dlt%sense (:,of) = 0.0_preal
344 dlt%chi0 (:,of) = 0.0_preal
345 dlt%gamma0(:,of) = 0.0_preal
350 end subroutine plastic_kinehardening_deltastate
356 module subroutine plastic_kinehardening_results(instance,group)
358 integer,
intent(in) :: instance
359 character(len=*),
intent(in) :: group
363 associate(prm => param(instance), stt => state(instance))
364 outputsloop:
do o = 1,
size(prm%output)
365 select case(trim(prm%output(o)))
367 if(prm%sum_N_sl>0)
call results_writedataset(group,stt%crss,
'xi_sl', &
368 'resistance against plastic slip',
'Pa')
370 if(prm%sum_N_sl>0)
call results_writedataset(group,stt%crss_back,
'tau_back', &
371 'back stress against plastic slip',
'Pa')
373 if(prm%sum_N_sl>0)
call results_writedataset(group,stt%sense,
'sense_of_shear', &
376 if(prm%sum_N_sl>0)
call results_writedataset(group,stt%chi0,
'chi0', &
379 if(prm%sum_N_sl>0)
call results_writedataset(group,stt%gamma0,
'gamma0', &
381 case (
'accumulatedshear')
382 if(prm%sum_N_sl>0)
call results_writedataset(group,stt%accshear,
'gamma_sl', &
388 end subroutine plastic_kinehardening_results
398 pure subroutine kinetics(Mp,instance,of, &
399 gdot_pos,gdot_neg,dgdot_dtau_pos,dgdot_dtau_neg)
401 real(pReal),
dimension(3,3),
intent(in) :: &
403 integer,
intent(in) :: &
407 real(pReal),
intent(out),
dimension(param(instance)%sum_N_sl) :: &
410 real(pReal),
intent(out),
optional,
dimension(param(instance)%sum_N_sl) :: &
414 real(pReal),
dimension(param(instance)%sum_N_sl) :: &
419 associate(prm => param(instance), stt => state(instance))
421 do i = 1, prm%sum_N_sl
422 tau_pos(i) = math_tensordot(mp,prm%nonSchmid_pos(1:3,1:3,i)) - stt%crss_back(i,of)
423 tau_neg(i) = merge(math_tensordot(mp,prm%nonSchmid_neg(1:3,1:3,i)) - stt%crss_back(i,of), &
424 0.0_preal, prm%nonSchmidActive)
427 where(dneq0(tau_pos))
428 gdot_pos = prm%gdot0 * merge(0.5_preal,1.0_preal, prm%nonSchmidActive) &
429 * sign(abs(tau_pos/stt%crss(:,of))**prm%n, tau_pos)
434 where(dneq0(tau_neg))
435 gdot_neg = prm%gdot0 * 0.5_preal &
436 * sign(abs(tau_neg/stt%crss(:,of))**prm%n, tau_neg)
441 if (
present(dgdot_dtau_pos))
then
442 where(dneq0(gdot_pos))
443 dgdot_dtau_pos = gdot_pos*prm%n/tau_pos
445 dgdot_dtau_pos = 0.0_preal
448 if (
present(dgdot_dtau_neg))
then
449 where(dneq0(gdot_neg))
450 dgdot_dtau_neg = gdot_neg*prm%n/tau_neg
452 dgdot_dtau_neg = 0.0_preal
457 end subroutine kinetics
459 end submodule plastic_kinehardening