forked from PolymerGuy/MaterialModel
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathumat.f
executable file
·409 lines (389 loc) · 19.6 KB
/
umat.f
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
SUBROUTINE VUMAT(
! READ ONLY - DO NOT MODIFY
. NBLOCK, NDIR, NSHR, NSTATEV, NFIELDV, NPROPS, LANNEAL, STEPTIME,
. TOTALTIME, DT, CMNAME, COORDMP, CHARLENGTH, PROPS, DENSITY,
. STRAININC, RELSPININC, TEMPOLD, STRETCHOLD, DEFGRADOLD, FIELDOLD,
. STRESSOLD, STATEOLD, ENERINTERNOLD, ENERINELASOLD, TEMPNEW,
. STRETCHNEW, DEFGRADNEW, FIELDNEW,
! WRITE ONLY - DO NOT READ
. STRESSNEW, STATENEW, ENERINTERNNEW, ENERINELASNEW)
!-----------------------------------------------------------------------
! ABAQUS implicit variable declaration included in VABA_PARAM.INC
! states the following:
! a to h are REAL*8 variables
! o to z are REAL*8 variables
! i to n are integer variables
!-----------------------------------------------------------------------
INCLUDE 'VABA_PARAM.INC'
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
! ABAQUS variables
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
DIMENSION PROPS(NPROPS), DENSITY(NBLOCK), COORDMP(NBLOCK,*),
. CHARLENGTH(*), STRAININC(NBLOCK,NDIR+NSHR), RELSPININC(*),
. TEMPOLD(*), STRETCHOLD(*), FIELDOLD(*),
. STRESSOLD(NBLOCK,NDIR+NSHR), STATEOLD(NBLOCK,NSTATEV),
. ENERINTERNOLD(NBLOCK), ENERINELASOLD(NBLOCK), TEMPNEW(*),
. STRETCHNEW(*), DEFGRADNEW(*), FIELDNEW(*),
. STRESSNEW(NBLOCK,NDIR+NSHR), STATENEW(NBLOCK,NSTATEV),
. ENERINTERNNEW(NBLOCK), ENERINELASNEW(NBLOCK)
C
CHARACTER*80 CMNAME
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
! Internal UMAT variables
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
DIMENSION STRESS(NDIR+NSHR) ! Stress tensor inside UMAT
DIMENSION DFDS(NDIR+NSHR) ! Derivative of the yield function
DIMENSION DFDS_SQRTJ2(NDIR+NSHR) ! Derivative of the yield function
DIMENSION STRESSK(NDIR+NSHR)
DIMENSION DSTRESSK(NDIR+NSHR)
DIMENSION EP(NDIR+NSHR)
DIMENSION DEP(NDIR+NSHR)
DIMENSION DEFGRADOLD(nblock,ndir+nshr+nshr)
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
! Internal UMAT variables could be declared as follow
! but it is not required due to the VABA_PARAM.INC file
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
! Material parameters
REAL*8 YOUNG ! Young's modulus
REAL*8 POISS ! Poisson's ratio
REAL*8 C11,C12,C44,CIIKK ! Elasticity matrix C components
REAL*8 SIGMA0 ! Initial yield stress
REAL*8 ET ! Tangent modulus
REAL*8 ALPHA ! Constant used in Drucker Prager model
! Visco plastic parameters
REAL*8 P0DOT ! Reference strain rate
REAL*8 S ! Visco plastic variable S
! Internal variables
REAL*8 P ! Equivalent plastic strain
! Plasticity variables
REAL*8 DLAMBDA ! Plastic multiplier
REAL*8 DDLAMBDA ! Increment in plastic multiplier
REAL*8 P0 ! Initial accumulated plastic strain
REAL*8 DPDT ! Plastic strain rate
! Visco plastic variables
REAL*8 SV ! Viscous back stress
REAL*8 DSVDP ! Partial derivative of viscous back stress
! Yield function variables
REAL*8 F ! Flow function
REAL*8 SQRTJ2 ! Equivalent von mises stress
REAL*8 I1 ! First stress invariant
REAL*8 PHI ! Equivalent stress
REAL*8 SIGMAY ! Yield stress
REAL*8 DFDP ! Gradient of "flow" function
! Computational variables
REAL*8 RESNOR ! Convergence criterion
REAL*8 TOL ! Tolerance for the RMAP algorithm
INTEGER MXITER ! Maximum number of iteration for the RMAP
INTEGER ITER ! Number of iteration for the RMAP
REAL*8 tst ! Variable for checking sign of function
! State variables
REAL*8 STRESSK ! Stress state at iteration k
! Constants
PARAMETER(three=3.0d0,two=2.0d0,one=1.0d0,zero=0.0d0)
ALPHA = 0.5d0
!-----------------------------------------------------------------------
! Read parameters from ABAQUS material card
!-----------------------------------------------------------------------
YOUNG = PROPS(1)
POISS = PROPS(2)
SIGMA0 = PROPS(3)
ET = PROPS(4)
S = PROPS(5)!*0.0d0
P0DOT = PROPS(6)
TOL = PROPS(7)
MXITER = PROPS(8)
!-----------------------------------------------------------------------
! Compute elasticity matrix
!-----------------------------------------------------------------------
C11 = YOUNG*(1.0d0-POISS)/((1.0d0+POISS)*(1.0d0-2.0d0*POISS)) !la+2nu
C12 = POISS*C11/(1.0d0-POISS) !Lame lambda
C44 = 0.5d0*YOUNG/(1.0d0+POISS) !Lame nu
CIIKK = three*YOUNG/(one-two*POISS)
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
! Loop over integration points
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
DO i=1,NBLOCK
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
! If time = 0 then pure elastic computation
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
IF(TOTALTIME.eq.0.0)THEN
STRESSNEW(i,1) = STRESSOLD(i,1)+C11*STRAININC(i,1)
. +C12*STRAININC(i,2)+C12*STRAININC(i,3)
STRESSNEW(i,2) = STRESSOLD(i,2)+C12*STRAININC(i,1)
. +C11*STRAININC(i,2)+C12*STRAININC(i,3)
STRESSNEW(i,3) = STRESSOLD(i,3)+C12*STRAININC(i,1)
. +C12*STRAININC(i,2)+C11*STRAININC(i,3)
STRESSNEW(i,4) = STRESSOLD(i,4)+C44*STRAININC(i,4)*2.0d0
STRESSNEW(i,5) = STRESSOLD(i,5)+C44*STRAININC(i,5)*2.0d0
STRESSNEW(i,6) = STRESSOLD(i,6)+C44*STRAININC(i,6)*2.0d0
ELSE
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
! Elastic-predictor-corrector scheme
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
!-----------------------------------------------------------------------
! Elastic prediction
!-----------------------------------------------------------------------
STRESS(1) = STRESSOLD(i,1)+C11*STRAININC(i,1)
. +C12*STRAININC(i,2)
. +C12*STRAININC(i,3)
STRESS(2) = STRESSOLD(i,2)+C12*STRAININC(i,1)
. +C11*STRAININC(i,2)
. +C12*STRAININC(i,3)
STRESS(3) = STRESSOLD(i,3)+C12*STRAININC(i,1)
. +C12*STRAININC(i,2)
. +C11*STRAININC(i,3)-DLAMBDA*C44*DFDS(6)
STRESS(4) = STRESSOLD(i,4)+C44*STRAININC(i,4)*2.0d0
STRESS(5) = STRESSOLD(i,5)+C44*STRAININC(i,5)*2.0d0
STRESS(6) = STRESSOLD(i,6)+C44*STRAININC(i,6)*2.0d0
STRESSK(1) = STRESS(1)
STRESSK(2) = STRESS(2)
STRESSK(3) = STRESS(3)
STRESSK(4) = STRESS(4)
STRESSK(5) = STRESS(5)
STRESSK(6) = STRESS(6)
!-----------------------------------------------------------------------
! Equivalent stress
!-----------------------------------------------------------------------
SQRTJ2 = sqrt(STRESSK(1)*STRESSK(1)
. +STRESSK(2)*STRESSK(2)
. +STRESSK(3)*STRESSK(3)
. -STRESSK(1)*STRESSK(2)
. -STRESSK(2)*STRESSK(3)
. -STRESSK(3)*STRESSK(1)
. +3.0d0*(STRESSK(4)*STRESSK(4)
. +STRESSK(5)*STRESSK(5)
. +STRESSK(6)*STRESSK(6)))
I1 = (STRESSK(1)+STRESSK(2)+STRESSK(3))
PHI = (SQRTJ2 + ALPHA*I1)/(one+ALPHA)
!-----------------------------------------------------------------------
! Equivalent plastic strain from previous time step
!-----------------------------------------------------------------------
P = STATEOLD(i,1)
P0 = STATEOLD(i,1)
DEP(1) = 0.0d0
DEP(2) = 0.0d0
DEP(3) = 0.0d0
DEP(4) = 0.0d0
DEP(5) = 0.0d0
DEP(6) = 0.0d0
!-----------------------------------------------------------------------
! Initialize the plastic multiplier
!-----------------------------------------------------------------------
DLAMBDA = 0.00000d0
!-----------------------------------------------------------------------
! Update yield stress
!-----------------------------------------------------------------------
SIGMAY = SIGMA0+ET*P
!-----------------------------------------------------------------------
! Compute yield function
!-----------------------------------------------------------------------
F = PHI-SIGMAY
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
! Check for plasticity
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
!-----------------------------------------------------------------------
! Plastic flow
!-----------------------------------------------------------------------
IF(F.GT.0.0d0)THEN ! Plastic flow
!-----------------------------------------------------------------------
! Determine viscous back stress
!-----------------------------------------------------------------------
SV = S*log(one+DLAMBDA/(DT*P0DOT))
DSVDP = S/(DT*P0DOT+DLAMBDA)
!-----------------------------------------------------------------------
! Compute augmented yield function
!-----------------------------------------------------------------------
F = PHI-SIGMAY - SV
!-----------------------------------------------------------------------
! Compute partial derivative of yield function
!-----------------------------------------------------------------------
DFDP = ET +(3.0d0*C44+(alpha**two)*CIIKK)/
. ((one+ALPHA)**two) +DSVDP
!-----------------------------------------------------------------------
! Initialize local iterations
!-----------------------------------------------------------------------
DO ITER=1,MXITER
!-----------------------------------------------------------------------
! Compute increment in plastic multiplier
!-----------------------------------------------------------------------
DDLAMBDA = F/DFDP
!-----------------------------------------------------------------------
! Update plastic multiplier and P
!-----------------------------------------------------------------------
DLAMBDA = DLAMBDA+DDLAMBDA
P = P0 + DLAMBDA
!-----------------------------------------------------------------------
! Update viscous back stress
!-----------------------------------------------------------------------
tst = one+DLAMBDA/(DT*P0DOT)
IF(tst.GT.zero)THEN
SV = S*log(tst)
DSVDP = (S/(tst))/(DT*P0DOT)
ELSE
SV = S*log(-tst)
DSVDP = -(S/(-tst))/(DT*P0DOT)
ENDIF
!-----------------------------------------------------------------------
! Compute the derivative of the yield function
!-----------------------------------------------------------------------
IF(PHI.eq.0)THEN
DENOM = 1.0d0
ELSE
DENOM = SQRTJ2
ENDIF
DFDS_SQRTJ2(1) = (STRESSK(1)-0.5d0*(STRESSK(2)+
. STRESSK(3)))/DENOM
DFDS_SQRTJ2(2) = (STRESSK(2)-0.5d0*(STRESSK(3)+
. STRESSK(1)))/DENOM
DFDS_SQRTJ2(3) = (STRESSK(3)-0.5d0*(STRESSK(1)+
. STRESSK(2)))/DENOM
DFDS_SQRTJ2(4) = 3.0d0*STRESSK(4)/(two*DENOM)
DFDS_SQRTJ2(5) = 3.0d0*STRESSK(5)/(two*DENOM)
DFDS_SQRTJ2(6) = 3.0d0*STRESSK(6)/(two*DENOM)
DFDS(1) = (DFDS_SQRTJ2(1)+ALPHA)/(one+ALPHA)
DFDS(2) = (DFDS_SQRTJ2(2)+ALPHA)/(one+ALPHA)
DFDS(3) = (DFDS_SQRTJ2(3)+ALPHA)/(one+ALPHA)
DFDS(4) = (DFDS_SQRTJ2(4))/(one+ALPHA)
DFDS(5) = (DFDS_SQRTJ2(5))/(one+ALPHA)
DFDS(6) = (DFDS_SQRTJ2(6))/(one+ALPHA)
!-----------------------------------------------------------------------
! Determine increment in plastic strain components
!-----------------------------------------------------------------------
DEP(1)= DLAMBDA*DFDS(1)
DEP(2)= DLAMBDA*DFDS(2)
DEP(3)= DLAMBDA*DFDS(3)
DEP(4)= DLAMBDA*DFDS(4)
DEP(5)= DLAMBDA*DFDS(5)
DEP(6)= DLAMBDA*DFDS(6)
!-----------------------------------------------------------------------
! Determine increment in stress components
!-----------------------------------------------------------------------
DSTRESSK(1) = C11*DEP(1)+C12*DEP(2)
. +C12*DEP(3)
DSTRESSK(2) = C12*DEP(1)+C11*DEP(2)
. +C12*DEP(3)
DSTRESSK(3) = C12*DEP(1)+C12*DEP(2)
. +C11*DEP(3)
DSTRESSK(4) = C44*two*DEP(4)
DSTRESSK(5) = C44*two*DEP(5)
DSTRESSK(6) = C44*two*DEP(6)
!-----------------------------------------------------------------------
! Update stress state
!-----------------------------------------------------------------------
STRESSK(1)= STRESS(1)-DSTRESSK(1)
STRESSK(2)= STRESS(2)-DSTRESSK(2)
STRESSK(3)= STRESS(3)-DSTRESSK(3)
STRESSK(4)= STRESS(4)-DSTRESSK(4)
STRESSK(5)= STRESS(5)-DSTRESSK(5)
STRESSK(6)= STRESS(6)-DSTRESSK(6)
!-----------------------------------------------------------------------
! Update yield stress
!-----------------------------------------------------------------------
SIGMAY = SIGMA0+ET*P
!-----------------------------------------------------------------------
! Equivalent stress
!-----------------------------------------------------------------------
SQRTJ2 = sqrt(STRESSK(1)*STRESSK(1)
. +STRESSK(2)*STRESSK(2)
. +STRESSK(3)*STRESSK(3)
. -STRESSK(1)*STRESSK(2)
. -STRESSK(2)*STRESSK(3)
. -STRESSK(3)*STRESSK(1)
. +3.0d0*(STRESSK(4)*STRESSK(4)
. +STRESSK(5)*STRESSK(5)
. +STRESSK(6)*STRESSK(6)))
I1 = STRESSK(1)+STRESSK(2)+STRESSK(3)
PHI = (SQRTJ2 + ALPHA*I1)/(one+ALPHA)
!-----------------------------------------------------------------------
! Compute augmented yield function
!-----------------------------------------------------------------------
F = PHI-SIGMAY - SV
!-----------------------------------------------------------------------
! Compute partial derivative of yield function
!-----------------------------------------------------------------------
DFDP = ET +(3.0d0*C44+(alpha**two)*CIIKK)/
. ((one+ALPHA)**two) +DSVDP
!-----------------------------------------------------------------------
! Compute convergence criterion
!-----------------------------------------------------------------------
RESNOR = ABS((F)/SIGMAY)
!-----------------------------------------------------------------------
! Check for convergence
!-----------------------------------------------------------------------
IF(RESNOR.LE.TOL.AND.DLAMBDA.GE.zero)THEN
PRINT *, 'Converged in', ITER
!-----------------------------------------------------------------------
! Update the stress tensor
!-----------------------------------------------------------------------
STRESSNEW(i,1) = STRESSK(1)
STRESSNEW(i,2) = STRESSK(2)
STRESSNEW(i,3) = STRESSK(3)
STRESSNEW(i,4) = STRESSK(4)
STRESSNEW(i,5) = STRESSK(5)
STRESSNEW(i,6) = STRESSK(6)
!-----------------------------------------------------------------------
! Update the history variables
!-----------------------------------------------------------------------
STATENEW(i,1) = P
STATENEW(i,2) = PHI-3.0*C44*DLAMBDA
STATENEW(i,3) = F
STATENEW(i,4) = SIGMAY
STATENEW(i,5) = DGAMA
STATENEW(i,6) = ITER
STATENEW(i,7) = STATEOLD(i,7) + DEP(1)
STATENEW(i,8) = STATEOLD(i,8) + DEP(2)
STATENEW(i,9) = STATEOLD(i,9) + DEP(3)
STATENEW(i,10) =STATEOLD(i,10) + DEP(4)
STATENEW(i,11) =STATEOLD(i,11) + DEP(5)
STATENEW(i,12) =STATEOLD(i,12) + DEP(6)
GOTO 90
ELSE ! RMAP has not converged yet
IF(ITER.eq.MXITER)THEN
write(*,*) 'RMAP has not converged'
write(*,*) 'Integration point',i
write(*,*) 'Convergence',RESNOR
write(*,*) 'dlambda',DLAMBDA,P
STOP
ENDIF
ENDIF
ENDDO
!-----------------------------------------------------------------------
! Elastic point
!-----------------------------------------------------------------------
ELSE
!-----------------------------------------------------------------------
! Update the stress tensor
!-----------------------------------------------------------------------
STRESSNEW(i,1) = STRESS(1)
STRESSNEW(i,2) = STRESS(2)
STRESSNEW(i,3) = STRESS(3)
STRESSNEW(i,4) = STRESS(4)
STRESSNEW(i,5) = STRESS(5)
STRESSNEW(i,6) = STRESS(6)
!-----------------------------------------------------------------------
! Update the history variables
!-----------------------------------------------------------------------
STATENEW(i,1) = P
STATENEW(i,2) = PHI
STATENEW(i,3) = F
STATENEW(i,4) = SIGMAY
STATENEW(i,5) = 0.0
STATENEW(i,6) = 0
STATENEW(i,7) = STATEOLD(i,7)
STATENEW(i,8) = STATEOLD(i,8)
STATENEW(i,9) = STATEOLD(i,9)
STATENEW(i,10) =STATEOLD(i,10)
STATENEW(i,11) =STATEOLD(i,11)
STATENEW(i,12) =STATEOLD(i,12)
ENDIF
!-----------------------------------------------------------------------
! End of loop over integration points
!-----------------------------------------------------------------------
ENDIF
90 CONTINUE
ENDDO
!-----------------------------------------------------------------------
! END SUBROUTINE
!-----------------------------------------------------------------------
RETURN
END