42
43
44
45
46
47
48 USE fail_param_mod
54
55
56
57#include "implicit_f.inc"
58
59
60
61#include "units_c.inc"
62
63
64
65 INTEGER ,INTENT(IN) :: FAIL_ID
66 INTEGER ,INTENT(IN) :: MAT_ID
67 INTEGER ,INTENT(IN) :: IRUPT
68 CHARACTER(LEN=NCHARTITLE) ,INTENT(IN) :: TITR
69 TYPE(UNIT_TYPE_) ,INTENT(IN) :: UNITAB
70 TYPE(SUBMODEL_DATA),INTENT(IN) :: LSUBMODEL(*)
71 TYPE(FAIL_PARAM_) ,INTENT(INOUT) :: FAIL
72
73
74
75 INTEGER :: NANGLE,I,J,K,INFO,REG_FUNC,MFLAG,SFLAG,RATE_FUNC,NFUNC,NUPARAM,NUVAR
76 INTEGER :: IPIV2(2),IPIV3(3)
77 INTEGER ,PARAMETER :: NSIZE = 2
78 INTEGER ,DIMENSION(NSIZE) :: IFUNC
79 my_real :: pthk,ref_siz,ref_siz_unit,epsd0,cjc,rate_scale,ref_rate_unit,
80 . r1,r2,r4,r5,c5,c5_min,theta_myreal
81 my_real,
DIMENSION(:),
ALLOCATABLE :: c1,c2,c3,c4,inst
82 DOUBLE PRECISION A_1(2,2),B_1(2),A_2(3,3),B_2(3),
83 . TRIAX_1_LIN,TRIAX_2_LIN,TRIAX_3_LIN,
84 . TRIAX_4_LIN,TRIAX_5_LIN,TRIAX_1_QUAD,
85 . TRIAX_2_QUAD,TRIAX_3_QUAD,TRIAX_4_QUAD,
86 . TRIAX_5_QUAD,COS2(10,10),XMIN,YMIN
87 DOUBLE PRECISION, DIMENSION(:), ALLOCATABLE :: THETA,THETA_RAD,Q_X11,Q_X12,Q_X13,
88 . Q_X21,Q_X22,Q_X23,Q_INST
89 DOUBLE PRECISION, DIMENSION(:,:), ALLOCATABLE :: X_1,X_2,AMAT,BVEC
90 INTEGER, DIMENSION(:), ALLOCATABLE :: IPIV
91 LOGICAL :: IS_AVAILABLE,IS_ENCRYPTED
92 DATA triax_1_lin, triax_2_lin, triax_3_lin, triax_4_lin,
93 . triax_5_lin
94 . / -0.33333333, 0.0, 0.33333333, 0.577350269, 0.66666667 /
95 DATA triax_1_quad, triax_2_quad, triax_3_quad,
96 . triax_4_quad, triax_5_quad
97 . / 0.11111111, 0.0, 0.11111111, 0.33333333, 0.44444444 /
98
99 DATA cos2/
100 1 1. ,0. ,0. ,0. ,0. ,0. ,0. ,0. ,0. ,0. ,
101 2 0. ,1. ,0. ,0. ,0. ,0. ,0. ,0. ,0. ,0. ,
102 3 -1. ,0. ,2. ,0. ,0. ,0. ,0. ,0. ,0. ,0. ,
103 4 0. ,-3. ,0. ,4. ,0. ,0. ,0. ,0. ,0. ,0. ,
104 5 1. ,0. ,-8. ,0. ,8. ,0. ,0. ,0. ,0. ,0. ,
105 6 0. ,5. ,0. ,-20. ,0. ,16. ,0. ,0. ,0. ,0. ,
106 7 -1. ,0. ,18. ,0. ,-48. ,0. ,32. ,0. ,0. ,0. ,
107 8 0. ,-7. ,0. ,56. ,0. ,-112.,0. ,64. ,0. ,0. ,
108 9 1. ,0. ,-32. ,0. ,160. ,0. ,-256.,0. ,128. ,0. ,
109 a 0. ,9. ,0. ,-120.,0. ,432. ,0. ,-576 ,0. ,256. /
110
111 is_encrypted = .false.
112 is_available = .false.
113
114
115
116
118
119
120
121
122
123
124
125 CALL hm_get_floatv (
'Pthk' ,pthk ,is_available,lsubmodel,unitab)
126 CALL hm_get_intv (
'MAT_MFLAG' ,mflag ,is_available,lsubmodel)
127 CALL hm_get_intv (
'MAT_SFLAG' ,sflag ,is_available,lsubmodel)
128 CALL hm_get_intv (
'MAT_refanglemax',nangle ,is_available,lsubmodel)
129
130 IF (nangle > 10) THEN
131 CALL ancmsg(msgid=2015,msgtype=msgerror,
132 . anmode=aninfo_blind_1,i1=mat_id,c1=titr)
133 ENDIF
134 CALL hm_get_intv (
'fct_IDel' ,reg_func ,is_available,lsubmodel)
135 CALL hm_get_floatv (
'EI_ref' ,ref_siz ,is_available,lsubmodel,unitab)
136
137 IF (pthk == zero) pthk = one - em06
138 pthk =
min(pthk, one)
139 pthk =
max(pthk,-one)
140 IF (sflag == 0) sflag = 2
141
142 IF ((ref_siz == zero).AND.(reg_func > 0)) THEN
144 ref_siz = one*ref_siz_unit
145 ENDIF
146
147
148
149
150 CALL hm_get_floatv (
'MAT_C5' ,c5 ,is_available,lsubmodel,unitab)
151 CALL hm_get_floatv (
'MAT_EPSD0' ,epsd0 ,is_available,lsubmodel,unitab)
152 CALL hm_get_floatv (
'MAT_CJC' ,cjc ,is_available,lsubmodel,unitab)
153 CALL hm_get_intv (
'fct_IDrate',rate_func ,is_available,lsubmodel)
154 CALL hm_get_floatv (
'RATE_scale',rate_scale ,is_available,lsubmodel,unitab)
155
156 IF ((rate_scale == zero).AND.(rate_func > 0)) THEN
157 CALL hm_get_floatv_dim(
'RATE_scale' ,ref_rate_unit ,is_available, lsubmodel, unitab)
158 rate_scale = ref_rate_unit*one
159 ENDIF
160 IF (rate_func > 0) THEN
161 cjc = zero
162 epsd0 = zero
163 ELSE
164 rate_scale = zero
165 ENDIF
166 IF (cjc == zero .OR. epsd0 == zero) THEN
167 cjc = zero
168 epsd0 = zero
169 ENDIF
170
171
172
173
174 IF (.NOT.ALLOCATED(c1)) ALLOCATE(c1(nangle))
175 IF (.NOT.ALLOCATED(c2)) ALLOCATE(c2(nangle))
176 IF (.NOT.ALLOCATED(c3)) ALLOCATE(c3(nangle))
177 IF (.NOT.ALLOCATED(c4)) ALLOCATE(c4(nangle))
178 IF (sflag == 3) THEN
179 IF (.NOT.ALLOCATED(inst)) ALLOCATE(inst(nangle))
180 inst = zero
181 ENDIF
182 c5_min = infinity
183
184 IF (mflag == 0) THEN
185
186 DO j = 1, nangle
191
192 IF (c3(j) == zero) c3(j) = 0.6d0
193 IF (c1(j) == zero .AND. c2(j) == zero .AND. c4(j) == zero .AND. c5 == zero) THEN
194 c1(j) = 3.5d0*c3(j)
195 c2(j) = 1.6d0*c3(j)
196 c4(j) = 0.6d0*c3(j)
197 c5_min =
min(c5_min,1.5d0*c3(j))
198 ENDIF
199
200 IF (sflag == 3) THEN
202 IF (inst(j) <= zero) THEN
203 CALL ancmsg(msgid=2016,msgtype=msgwarning,
204 . anmode=aninfo_blind_1,i1=mat_id,c1=titr)
205 sflag = 2
206 ELSEIF (inst(j) >= c4(j)) THEN
207 CALL ancmsg(msgid=2017,msgtype=msgwarning,
208 . anmode=aninfo_blind_1,i1=mat_id,c1=titr)
209 sflag = 2
210 ENDIF
211 ENDIF
212 ENDDO
213 ELSE
214
215 IF (mflag == 1) THEN
216 r1 = 3.5d0
217 r2 = 1.6d0
218 r4 = 0.6d0
219 r5 = 1.5d0
220 ELSEIF (mflag == 2) THEN
221 r1 = 4.3d0
222 r2 = 1.4d0
223 r4 = 0.6d0
224 r5 = 1.6d0
225 ELSEIF (mflag == 3) THEN
226 r1 = 5.2d0
227 r2 = 3.1d0
228 r4 = 0.8d0
229 r5 = 3.5d0
230 ELSEIF (mflag == 4) THEN
231 r1 = 5.0d0
232 r2 = 1.0d0
233 r4 = 0.4d0
234 r5 = 0.8d0
235 ELSEIF (mflag == 5) THEN
236 r1 = 7.8d0
237 r2 = 3.5d0
238 r4 = 0.6d0
239 r5 = 2.8d0
240 ELSEIF (mflag == 6) THEN
241 r1 = 3.6d0
242 r2 = 0.6d0
243 r4 = 0.5d0
244 r5 = 0.6d0
245 ELSEIF (mflag == 7) THEN
246 r1 = 10.0d0
247 r2 = 2.7d0
248 r4 = 0.6d0
249 r5 = 0.7d0
250 ELSEIF (mflag == 99) THEN
251 CALL hm_get_floatv (
'MAT_R1' ,r1 ,is_available,lsubmodel,unitab)
252 CALL hm_get_floatv (
'MAT_R2' ,r2 ,is_available,lsubmodel,unitab)
253 CALL hm_get_floatv (
'MAT_R4' ,r4 ,is_available,lsubmodel,unitab)
254 CALL hm_get_floatv (
'MAT_R5' ,r5 ,is_available,lsubmodel,unitab)
255 ELSE
256 r1 = 3.5d0
257 r2 = 1.6d0
258 r4 = 0.6d0
259 r5 = 1.5d0
260 ENDIF
261
262 DO j = 1, nangle
264
265 IF (c3(j) == zero) THEN
266 IF (mflag == 1) THEN
267 c3(j) = 0.6d0
268 ELSEIF (mflag == 2) THEN
269 c3(j) = 0.5d0
270 ELSEIF (mflag == 3) THEN
271 c3(j) = 0.12d0
272 ELSEIF (mflag == 4) THEN
273 c3(j) = 0.3d0
274 ELSEIF (mflag == 5) THEN
275 c3(j) = 0.17d0
276 ELSEIF (mflag == 6) THEN
277 c3(j) = 0.1d0
278 ELSEIF (mflag == 7) THEN
279 c3(j) = 0.11d0
280 ENDIF
281 ENDIF
282
283 c1(j) = r1*c3(j)
284 c2(j) = r2*c3(j)
285 c4(j) = r4*c3(j)
286 c5_min =
min(c5_min,r5*c3(j))
287
288 IF (sflag == 3) THEN
290 IF (inst(j) <= zero) THEN
291 CALL ancmsg(msgid=2016,msgtype=msgwarning,
292 . anmode=aninfo_blind_1,i1=mat_id,c1=titr)
293 sflag = 2
294 ELSEIF (inst(j) >= c4(j)) THEN
295 CALL ancmsg(msgid=2017,msgtype=msgwarning,
296 . anmode=aninfo_blind_1,i1=mat_id,c1=titr)
297 sflag = 2
298 ENDIF
299 ENDIF
300 ENDDO
301 ENDIF
302
303 IF (c5 == zero) c5 = c5_min
304
305
306
307 IF (.NOT.ALLOCATED(x_1)) ALLOCATE(x_1(nangle,3))
308 IF (.NOT.ALLOCATED(x_2)) ALLOCATE(x_2(nangle,3))
309
310
311
312 DO j = 1,nangle
313
314
315
316 a_1(1,1) = triax_1_lin
317 a_1(1,2) = triax_1_quad
318 a_1(2,1) = triax_3_lin
319 a_1(2,2) = triax_3_quad
320 b_1(1) = c1(j) - c2(j)
321 b_1(2) = c3(j) - c2(j)
322
323
324#ifndef WITHOUT_LINALG
325 CALL dgesv(2, 1, a_1, 2, ipiv2, b_1, 2, info)
326#else
327 WRITE(6,*) "Error: Blas/Lapack required"
328#endif
329 x_1(j,1) = c2(j)
330 x_1(j,2:3) = b_1(1:2)
331
332
333
334 a_2(1,1) = one
335 a_2(1,2) = triax_3_lin
336 a_2(1,3) = triax_3_quad
337 a_2(2,1) = one
338 a_2(2,2) = triax_4_lin
339 a_2(2,3) = triax_4_quad
340 a_2(3,1) = one
341 a_2(3,2) = triax_5_lin
342 a_2(3,3) = triax_5_quad
343 b_2(1) = c3(j)
344 b_2(2) = c4(j)
345 b_2
346
347
348#ifndef WITHOUT_LINALG
349 CALL dgesv(3, 1, a_2, 3, ipiv3, b_2, 3, info)
350#endif
351 x_2(j,1:3) = b_2(1:3)
352
353 ENDDO
354
355
356
357
358
359 IF (.NOT.ALLOCATED(theta)) ALLOCATE(theta(nangle))
360 IF (.NOT.ALLOCATED(theta_rad)) ALLOCATE(theta_rad(nangle))
361
362
363 DO j = 1, nangle
364 IF (nangle > 1) THEN
365 theta(j) = (j-1)*(90.0d0/(nangle-1))
366 theta_rad(j) = theta(j)*(pi/180.0d0)
367 ELSE
368 theta(j) = zero
369 theta_rad(j) = zero
370 ENDIF
371
372
373 xmin = -x_1(j,2)/(two*x_1(j,3))
374 ymin = x_1(j,3)*(xmin**2) + x_1(j,2)*xmin + x_1(j,1)
375 IF (ymin < zero) THEN
376 theta_myreal = theta(j)
378 . msgtype=msgwarning,
379 . anmode=aninfo_blind_1,
380 . i1=mat_id,
381 . c1=titr,
382 . i2=j,
383 . r1=theta_myreal)
384 ENDIF
385
386
387 IF (sflag == 1) THEN
388 xmin = -x_2(j,2)/(two*x_2(j,3))
389 ymin = x_2(j,3)*(xmin**2) + x_2(j,2)*xmin + x_2(j,1)
390 IF (ymin < zero) THEN
391 theta_myreal = theta(j)
393 . msgtype=msgwarning,
394 . anmode=aninfo_blind_1,
395 . i1=mat_id,
396 . c1=titr,
397 . i2=j,
398 . r1=theta_myreal)
399 ENDIF
400 ENDIF
401
402 ENDDO
403
404
405 ALLOCATE (amat(nangle,nangle),ipiv(nangle))
406
407
408 DO j = 1,nangle
409 DO i = 1,nangle
410 amat(j,i) = zero
411 DO k = 1,i
412 amat(j,i) = amat(j,i) + cos2(k,i)*(cos(two*theta_rad(j)))**(k-1)
413 ENDDO
414 ENDDO
415 ENDDO
416
417
418 ALLOCATE(q_x11(nangle),q_x12
419 . q_x21(nangle),q_x22(nangle),q_x23(nangle))
420
421
422 q_x11(1:nangle) = zero
423 q_x12(1:nangle) = zero
424 q_x13(1:nangle) = zero
425 q_x21(1:nangle) = zero
426 q_x22(1:nangle) = zero
427 q_x23(1:nangle) = zero
428 IF (sflag == 3) THEN
429 ALLOCATE(q_inst(nangle))
430 q_inst(1:nangle) = zero
431 ENDIF
432
433
434 IF (sflag == 3) THEN
435 ALLOCATE (bvec(nangle,7))
436 ELSE
437 ALLOCATE (bvec(nangle,6))
438 ENDIF
439 bvec(1:nangle,1) = x_1(1:nangle,1)
440 bvec(1:nangle,2) = x_1(1:nangle,2)
441 bvec(1:nangle,3) = x_1(1:nangle,3)
442 bvec(1:nangle,4) = x_2(1:nangle,1)
443 bvec(1:nangle,5) = x_2(1:nangle,2)
444 bvec(1:nangle,6) = x_2(1:nangle,3)
445 IF (sflag == 3) bvec(1:nangle,7) = inst(1:nangle)
446
447
448 ipiv(1:nangle) = 0
449
450
451#ifndef WITHOUT_LINALG
452 IF (sflag == 3) THEN
453 CALL dgesv(nangle, 7, amat, nangle, ipiv, bvec, nangle, info)
454 ELSE
455 CALL dgesv(nangle, 6, amat, nangle, ipiv, bvec, nangle, info)
456 ENDIF
457#else
458 WRITE(6,*) "Error: Blas/Lapack required"
459#endif
460
461
462 q_x11(1:nangle) = bvec(1:nangle,1)
463 q_x12(1:nangle) = bvec(1:nangle,2)
464 q_x13(1:nangle) = bvec(1:nangle,3)
465 q_x21(1:nangle) = bvec(1:nangle,4)
466 q_x22(1:nangle) = bvec(1:nangle,5)
467 q_x23(1:nangle) = bvec(1:nangle,6)
468 IF (sflag == 3) q_inst(1:nangle) = bvec(1:nangle,7)
469
470
471
472 nuparam = 7
473 IF (sflag == 3) THEN
474 nuparam = nuparam + 7*nangle
475 ELSE
476 nuparam = nuparam + 6*nangle
477 ENDIF
478
479 nfunc = 0
480 IF (rate_func /= 0) THEN
481 nfunc = nfunc + 1
482 ifunc(nfunc) = rate_func
483 ENDIF
484 IF (reg_func /= 0) THEN
485 nfunc = nfunc + 1
486 ifunc(nfunc) = reg_func
487 ENDIF
488
489 nuvar = 3
490
491
492
493 fail%KEYWORD = 'ORTH-BIQUAD'
494 fail%IRUPT = irupt
495 fail%FAIL_ID = fail_id
496 fail%NUPARAM = nuparam
497 fail%NIPARAM = 0
498 fail%NUVAR = nuvar
499 fail%NFUNC = nfunc
500 fail%NTABLE = 0
501 fail%NMOD = 0
502 fail%PTHK = pthk
503
504 ALLOCATE (fail%UPARAM(fail%NUPARAM))
505 ALLOCATE (fail%IPARAM(fail%NIPARAM))
506 ALLOCATE (fail%IFUNC (fail%NFUNC))
507 ALLOCATE (fail%TABLE (fail%NTABLE))
508
509 fail%IFUNC(1:nfunc) = ifunc(1:nfunc)
510
511 fail%UPARAM(1) = pthk
512 fail%UPARAM(2) = sflag
513 fail%UPARAM(3) = ref_siz
514 fail%UPARAM(4) = epsd0
515 fail%UPARAM(5) = cjc
516 fail%UPARAM(6) = rate_scale
517 fail%UPARAM(7) = nangle
518 IF (sflag == 3) THEN
519 DO j = 1,nangle
520 fail%UPARAM(8 + 7*(j-1)) = q_x11(j)
521 fail%UPARAM(9 + 7*(j-1)) = q_x12(j)
522 fail%UPARAM(10 + 7*(j-1)) = q_x13(j)
523 fail%UPARAM(11 + 7*(j-1)) = q_x21(j)
524 fail%UPARAM(12 + 7*(j-1)) = q_x22(j)
525 fail%UPARAM(13 + 7*(j-1)) = q_x23(j)
526 fail%UPARAM(14 + 7*(j-1)) = q_inst(j)
527 ENDDO
528 ELSE
529 DO j = 1,nangle
530 fail%UPARAM(8 + 6*(j-1)) = q_x11(j)
531 fail%UPARAM(9 + 6*(j-1)) = q_x12(j)
532 fail%UPARAM(10 + 6*(j-1)) = q_x13(j)
533 fail%UPARAM(11 + 6*(j-1)) = q_x21(j)
534 fail%UPARAM(12 + 6*(j-1)) = q_x22(j)
535 fail%UPARAM(13 + 6*(j-1)) = q_x23(j)
536 ENDDO
537 ENDIF
538
539
540
541 IF (is_encrypted) THEN
542 WRITE(iout,'(5X,A,//)')'CONFIDENTIAL DATA'
543 ELSE
544 WRITE(iout,1000)
545 IF (mflag /= 0) WRITE(iout, 1100) mflag
546 DO j = 1,nangle
547 WRITE(iout,1200) theta(j),c1(j),c2(j),c3(j),c4(j),c5,
548 & x_1(j,3),x_1(j,2),x_1(j,1),x_2(j,3),x_2(j,2),x_2(j,1)
549 IF (sflag == 3) WRITE(iout, 1900) inst(j)
550 ENDDO
551 WRITE(iout,1300) pthk,sflag
552 IF (reg_func > 0) WRITE(iout,1400) reg_func,ref_siz
553 IF (epsd0 > zero) THEN
554 WRITE(iout,1500) epsd0,cjc
555 ELSEIF (rate_func > 0) THEN
556 WRITE(iout,1600) rate_func,rate_scale
557 ENDIF
558 WRITE(iout,2000)
559 ENDIF
560
561
562
563 IF (ALLOCATED(c1)) DEALLOCATE(c1)
564 IF (ALLOCATED(c2)) DEALLOCATE(c2)
565 IF (ALLOCATED(c3)) DEALLOCATE(c3)
566 IF (ALLOCATED(c4)) DEALLOCATE(c4)
567 IF (ALLOCATED(inst)) DEALLOCATE(inst)
568 IF (ALLOCATED(x_1)) DEALLOCATE(x_1)
569 IF (ALLOCATED(x_2)) DEALLOCATE(x_2)
570 IF (ALLOCATED(theta)) DEALLOCATE(theta)
571 IF (ALLOCATED(theta_rad)) DEALLOCATE(theta_rad)
572 IF (ALLOCATED(q_x11)) DEALLOCATE(q_x11)
573 IF (ALLOCATED(q_x12)) DEALLOCATE(q_x12)
574 IF (ALLOCATED(q_x13)) DEALLOCATE(q_x13)
575 IF (ALLOCATED(q_x21)) DEALLOCATE(q_x21)
576 IF (ALLOCATED(q_x22)) DEALLOCATE(q_x22)
577 IF (ALLOCATED(q_x23)) DEALLOCATE(q_x23)
578 IF (ALLOCATED(q_inst)) DEALLOCATE(q_inst)
579 IF (ALLOCATED(amat)) DEALLOCATE(amat)
580 IF (ALLOCATED(ipiv)) DEALLOCATE(ipiv)
581
582 1000 FORMAT(
583 & 5x,' ------------------------------------------ ',/
584 & 5x,' FAILURE CRITERION : ORTHOTROPIC BIQUAD ',/,
585 & 5x,' ------------------------------------------ ',/)
586 1100 FORMAT(
587 & 5x,'MATERIAL PARAMETER SELECTOR . . . . . . . . . . .=',i10)
588 1200 FORMAT(
589 & 5x,'|| FAILURE STRAINS FOR ANGLE',f5.1,' DEG',/,
590 & 5x,' -------------------------------------------------',/,
591 & 5x,' SIMPLE COMPRESSION C1 . . . . . . . . . . . . .=',1pg20.13/
592 & 5x,' SHEAR C2 . . . . . . . . . . . . . . . . . . . .=',1pg20.13/
593 & 5x,' SIMPLE TENSION C3 . . . . . . . . . . . . . . .='
594 & 5x,' PLANE STRAIN C4 . . . . . . . . . . . . . . . .=',1pg20.13/
595 & 5x,' BIAXIAL TENSION C5 . . . . . . . . . . . . . . .=',1pg20.13/
596 & 5x,' ',/
597 & 5x,' LOW STRESS TRIAXIALITY PARABOLA PARAMETER A. . .='
598 & 5x,' LOW STRESS TRIAXIALITY PARABOLA PARAMETER B. . .=',1pg20.13/
599 & 5x,' LOW STRESS TRIAXIALITY PARABOLA PARAMETER C. . .=',1pg20.13/
600 & 5x,' ',/
601 & 5x,' high stress triaxiality parabola PARAMETER d . .=',1PG20.13/
602 & 5X,' high stress triaxiality parabola PARAMETER e . .=',1PG20.13/
603 & 5X,' high stress triaxiality parabola PARAMETER f . .=',1PG20.13/)
604 1300 FORMAT(
605 & 5X,'element deletion :',/,
606 & 5X,'shell element deletion PARAMETER pthickfail. . . .=',1PG20.13,/,
607 & 5X,' > 0.0 : fraction of failed thickness ',/,
608 & 5X,' < 0.0 : fraction of failed intg. points or layers',/,
609 & 5X,'s-flag . . . . . . . . . . . . . . . . . . . . . .=',I10,/)
610 1400 FORMAT(
611 & 5X,'element length regularization used:',/,
612 & 5X,'regularization
FUNCTION id . . . . . . . . . . . .=
',I10,/,
613 & 5X,'reference element length . . . . . . . . . . . . .=',1PG20.13,/)
614 1500 FORMAT(
615 & 5X,'johnson-cook strain-rate dependency:',/,
616 & 5X,'reference strain-rate . . . . . . . . . . . . . .=',1PG20.13,/,
617 & 5X,'johnson-cook parameter . . . . . . . . . . . . . .=',1PG20.13,/)
618 1600 FORMAT(
619 & 5X,'tabulated strain-rate dependency:
',/,
620 & 5X,'strain-rate dependency function
id . . . . . . . .=
',I10,/,
621 & 5X,'strain-rate scale factor . . . . . . . . . . . . .=',1PG20.13,/)
622 1900 FORMAT(
623 & 5X,' instability strain . . . . . . . . . . . . . . .=',1PG20.13,//)
624 2000 FORMAT(
625 & 5X,' --------------------------------------------------',//)
626
subroutine dgesv(n, nrhs, a, lda, ipiv, b, ldb, info)
DGESV computes the solution to system of linear equations A * X = B for GE matrices
subroutine hm_get_float_array_index(name, rval, index, is_available, lsubmodel, unitab)
subroutine hm_get_floatv(name, rval, is_available, lsubmodel, unitab)
subroutine hm_get_floatv_dim(name, dim_fac, is_available, lsubmodel, unitab)
subroutine hm_get_intv(name, ival, is_available, lsubmodel)
subroutine hm_option_is_encrypted(is_encrypted)
integer, parameter nchartitle
subroutine ancmsg(msgid, msgtype, anmode, i1, i2, i3, i4, i5, i6, i7, i8, i9, i10, i11, i12, i13, i14, i15, i16, i17, i18, i19, i20, r1, r2, r3, r4, r5, r6, r7, r8, r9, c1, c2, c3, c4, c5, c6, c7, c8, c9, prmode)
subroutine tabulated(iflag, nel, pm, off, eint, mu, espe, dvol, df, vnew, mat, psh, pnew, dpdm, dpde, npf, tf)