-
Notifications
You must be signed in to change notification settings - Fork 2
/
StressSolve.f90
878 lines (711 loc) · 27.6 KB
/
StressSolve.f90
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
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
INCLUDE 'IntegrateCrystalEqns.f90'
!***************************************
SUBROUTINE StressSolve( &
stress, ddsdde,strain, &
dstrain, statev, nstatv, &
time, dtime, kinc, kstep, &
pnewdt, argmin, argmax, &
geqvalues, &
gstress_n, &
gestran_n, &
gkappa_n, &
gmu, &
gmu_n, &
gdd_for_n, &
gdd_revp_n, &
gdd_revn_n, &
gdd_stressf_n, &
gdd_rev0_n, &
gdd_deb_n, &
noel, &
grainid)
!***************************************
USE PlaPar
USE NumType
USE FILEIO
USE SlipGeo
USE IterPar
USE NbData
IMPLICIT NONE
INTEGER :: nstatv ,kinc, kstep, iNG, noel
REAL(KIND=8) :: argmin, argmax, dtime, pnewdt, time(2)
REAL(KIND=8) :: stress(NVEC), ddsdde(NVEC, NVEC), PartStress(NVEC), &
strain(NVEC), dstrain(NVEC), statev(nstatv)
REAL(KIND=8) :: gstress_n (NVEC)
REAL(KIND=8) :: gstress (NVEC)
REAL(KIND=8) :: stress_n (NVEC)
REAL(KIND=8) :: estran (NVEC)
REAL(KIND=8) :: gestran (NVEC)
REAL(KIND=8) :: estran_n (NVEC)
REAL(KIND=8) :: gestran_n (NVEC)
REAL(KIND=8) :: geqvalues (NEQVA)
REAL(KIND=8) :: eqvalues (NEQVA)
REAL(KIND=8) :: gmu (NVEC)
REAL(KIND=8) :: gmu_n (NVEC)
REAL(KIND=8) :: mu (NVEC)
REAL(KIND=8) :: mu_n (NVEC)
REAL(KIND=8) :: gkappa_n (MaxNumSlip)
REAL(KIND=8) :: kappa (MaxNumSlip)
REAL(KIND=8) :: gkappa (MaxNumSlip)
REAL(KIND=8) :: kappa_n (MaxNumSlip)
REAL(KIND=8) :: gdd_for_n (MaxNumSlip)
REAL(KIND=8) :: gdd_revp_n (MaxNumSlip)
REAL(KIND=8) :: gdd_revn_n (MaxNumSlip)
REAL(KIND=8) :: gdd_stressf_n (MaxNumSlip)
REAL(KIND=8) :: gdd_rev0_n (MaxNumSlip)
REAL(KIND=8) :: gdd_deb_n (MaxNumSlip)
REAL(KIND=8) :: gdd_for (MaxNumSlip)
REAL(KIND=8) :: gdd_revp (MaxNumSlip)
REAL(KIND=8) :: gdd_revn (MaxNumSlip)
REAL(KIND=8) :: gdd_stressf (MaxNumSlip)
REAL(KIND=8) :: gdd_rev0 (MaxNumSlip)
REAL(KIND=8) :: gdd_deb (MaxNumSlip)
REAL(KIND=8) :: dd_for (MaxNumSlip)
REAL(KIND=8) :: dd_revp (MaxNumSlip)
REAL(KIND=8) :: dd_revn (MaxNumSlip)
REAL(KIND=8) :: dd_stressf (MaxNumSlip)
REAL(KIND=8) :: dd_rev0 (MaxNumSlip)
REAL(KIND=8) :: dd_deb (MaxNumSlip)
REAL(KIND=8) :: dd_for_n (MaxNumSlip)
REAL(KIND=8) :: dd_revp_n (MaxNumSlip)
REAL(KIND=8) :: dd_revn_n (MaxNumSlip)
REAL(KIND=8) :: dd_stressf_n (MaxNumSlip)
REAL(KIND=8) :: dd_rev0_n (MaxNumSlip)
REAL(KIND=8) :: dd_deb_n (MaxNumSlip)
REAL(KIND=8) :: ggamdot (MaxNumSlip)
REAL(KIND=8) :: gamdot (MaxNumSlip)
REAL(KIND=8) :: grss (MaxNumSlip)
REAL(KIND=8) :: rss (MaxNumSlip)
REAL(KIND=8) :: gmudot (NVEC)
REAL(KIND=8) :: mudot (NVEC)
REAL(KIND=8) :: dd_all
REAL(KIND=8) :: dd_max
! history maximum dislocation density
REAL(KIND=8) :: dd_his
REAL(KIND=8) :: dd_his_max
REAL(KIND=8) :: det_dd, det_dd_tmp
REAL(KIND=8) :: det_dd_max
REAL(KIND=8) :: ave_estran (NVEC)
REAL(KIND=8) :: det_dd_max_his
REAL(KIND=8) :: savg_ij(NSD, NSD)
REAL(KIND=8) :: cavg_ijkl(NVEC, NVEC)
REAL(KIND=8) :: c_ijkl(NVEC, NVEC)
REAL(KIND=8) :: s_ij(NSD, NSD)
INTEGER statusFlag
REAL(KIND=8) :: VecM(NSD,MaxNumSlip)
REAL(KIND=8) :: VecS(NSD,MaxNumSlip)
REAL(KIND=8) :: ZVec(NVEC,MaxNumSlip)
INTEGER :: grainid,is, IPart1, IPart2, ierr, iterState
CHARACTER(LEN=30) errinfo
!-----------------------------------------------------------------------
if (debug==1) then
write(FILE_E,*) '====***==== Begin ROMSolve ====***===='
Write(FILE_E,*) 'Dstrain and strain is:'
Write(FILE_E,'(6e18.8)') dstrain
Write(FILE_E,'(6e18.8)') strain
write(FILE_E,1000) time(1),dtime,kinc
endif
! fetch state variables from abaqus vector array
call RecoverStateVars(statev, nstatv, gstress_n, gestran_n, &
gkappa_n, ggamdot, grss, gmudot,geqvalues, gmu ,gmu_n, &
gdd_for_n, gdd_for, gdd_revp_n, gdd_revp, gdd_revn_n, &
gdd_revn,gdd_stressf_n, gdd_stressf, gdd_rev0_n, gdd_rev0, &
gdd_deb_n, gdd_deb, grainid)
! compute state
statusFlag = XTAL_CONVERGED
! Initialize the streess of each part and system Jacobian
savg_ij=pzero
cavg_ijkl=pzero
! Initialize some local arrays
stress=pzero
PartStress=pzero
estran=pzero
kappa=pzero
mu=pzero
c_ijkl=pzero
rss=pzero
mudot=pzero
gamdot=pzero
stress_n=pzero
kappa_n=pzero
mu_n=pzero
dd_for_n=pzero
dd_revp_n=pzero
dd_revn_n=pzero
dd_stressf_n=pzero
dd_rev0_n=pzero
dd_deb_n=pzero
dd_all=pzero
dd_max=pzero
! Initialize the history dd max and max of dd discrepancy
dd_his=statev(nstatv-10)! dd_his_max in each part
dd_his_max=statev(nstatv-8) ! dd_history_max in all parts
ave_estran(:) = statev((nstatv-6):(nstatv-1)) ! average plastic strain in each part
det_dd_max_his=statev(nstatv) ! dd_history_det_max_DD in all parts
! Fetch Crystal Variables At Part
if (debug==1) then
Write(FILE_E,*) 'initial estimate stress_n='
Write(FILE_E,*) gstress_n
end if
call FetchCrystalVariablesAtPart(gstress_n, gestran_n, gkappa_n, &
geqvalues,stress_n, estran_n, kappa_n, eqvalues, mu_n, &
gmu_n, dd_for_n, gdd_for_n, dd_revp_n, gdd_revp_n, &
dd_revn_n, gdd_revn_n, dd_stressf_n, gdd_stressf_n, &
dd_rev0_n, gdd_rev0_n, dd_deb_n, gdd_deb_n)
if (debug==1) then
write(FILE_E,*) 'initial estimate stress_n='
write(FILE_E,*) stress_n
end if
call IntegrateCrystalEqns2(s_ij, partstress, estran, kappa, mu, &
mu_n, eqvalues, gamdot, rss, mudot, stress_n, estran_n, &
kappa_n, dtime, iterP%maxIterstate, iterP%MaxIterNewt, ierr, &
iterState, iterP%tolerState, iterP%tolerNewt, c_ijkl, dstrain, kinc, &
kstep, argmin, argmax, SlipG%VecM0(:,:,grainid), SlipG%VecS0(:,:,grainid), SlipG%ZVEC0(:,:, grainid),dd_for_n, dd_for, &
dd_revp_n, dd_revp, dd_revn_n, dd_revn, dd_stressf_n, dd_stressf, &
dd_rev0_n, dd_rev0, dd_deb_n, dd_deb, time(2), dd_all, dd_max, &
noel, grainid)
if (ierr==XTAL_CONVERGED) then
errinfo='Converged'
elseif (ierr==XTAL_SING_JACOBIAN) then
errinfo= 'SING_JACOBIAN'
elseif (ierr==XTAL_MAX_ITERS_HIT) then
errinfo= 'MAX_ITERS_HIT'
else
errinfo= 'unknown error message, impossible!'
call RuntimeError(FILE_O, 'unknown error message, impossible!')
endif
if (iterprint==1) write(iter_O, '(3(I6,6x), 6x, A20)'), kstep, kinc, iterstate, errinfo
if (debug ==1) write(FILE_E,*) 'ierr= ', ierr
IF (ierr .ne. XTAL_CONVERGED) then
call WriteMessage(FILE_O, 'Resetting xtal quantities')
call RecoverStateVars(statev, nstatv, gstress_n, gestran_n, &
gkappa_n, ggamdot, grss, gmudot,geqvalues, gmu ,gmu_n, &
gdd_for_n, gdd_for, gdd_revp_n, gdd_revp, gdd_revn_n, &
gdd_revn,gdd_stressf_n, gdd_stressf, gdd_rev0_n, gdd_rev0, &
gdd_deb_n, gdd_deb, grainid)
call FetchCrystalVariablesAtPart(gstress_n, gestran_n, gkappa_n, &
geqvalues,stress_n, estran_n, kappa_n, eqvalues, mu_n, &
gmu_n, dd_for_n, gdd_for_n, dd_revp_n, gdd_revp_n, &
dd_revn_n, gdd_revn_n, dd_stressf_n, gdd_stressf_n, &
dd_rev0_n, gdd_rev0_n, dd_deb_n, gdd_deb_n)
call ResetCrystalQnts(partstress, estran, kappa, statev, eqvalues, &
gamdot, rss, mudot, mu, mu_n, s_ij, c_ijkl, stress_n, &
estran_n, kappa_n, dtime, SlipG%VecM0, SlipG%VecS0, SlipG%ZVEC0, dd_for_n,&
dd_for, dd_revp_n, dd_revp, dd_revn_n, dd_revn, &
dd_stressf_n, dd_stressf, dd_rev0_n, dd_rev0, dd_deb_n, &
dd_deb, argmin, argmax, grainid )
call SaveCrystalVariablesAtPart(partstress, estran, kappa, gamdot, &
rss, mudot, eqvalues, gstress, gestran, gkappa, gmu, mu, &
ggamdot, grss, gmudot, geqvalues, dd_for, gdd_for, dd_revp, &
gdd_revp, dd_revn, gdd_revn, dd_stressf, gdd_stressf, &
dd_rev0, gdd_rev0, dd_deb, gdd_deb)
call SaveStateVars(statev, nstatv, gstress, gestran, gkappa, &
ggamdot, grss, gmudot, gmu, geqvalues, gdd_for, gdd_revp, &
gdd_revn, gdd_stressf, gdd_rev0, gdd_deb, dd_all, dd_max, &
dd_his, dd_his_max, det_dd_max, ave_estran, det_dd_max_his, &
grainid)
write(FILE_O, *) ' ** Umat did not converged **'
write(FILE_O, *) ' ** re-scaling time step by 0.75 **'
pnewdt = 0.75
savg_ij=s_ij
cavg_ijkl=c_ijkl
CALL SaveStressModuli(stress, ddsdde, savg_ij, cavg_ijkl)
RETURN
ENDIF
call SaveCrystalVariablesAtPart(partstress, estran, kappa, gamdot, &
rss, mudot, eqvalues, gstress, gestran, gkappa, gmu, mu, &
ggamdot, grss, gmudot, geqvalues, dd_for, gdd_for, dd_revp, &
gdd_revp, dd_revn, gdd_revn, dd_stressf, gdd_stressf, &
dd_rev0, gdd_rev0, dd_deb, gdd_deb)
savg_ij=s_ij
cavg_ijkl=c_ijkl
!***********************************************************************
!-->damin: comment this part out, incompatible with direct simulation
! Output the max of dd discrepancy and local history max
! det_dd_tmp=0.d0
! if (dd_all .ge. dd_his) dd_his=dd_all
! DO iNG = 1, Elset_nb(grainid)%NNG
! det_dd_tmp(Elset_nb(grainid)%NBGID(iNG))=abs(dd_all(grianid)-dd_all(Elset_nb(grainid)%NBGID(iNG)))
! ENDDO
! det_dd=maxval(det_dd_tmp)
! det_dd_max=maxval(det_dd)
!! Output the max inside history dd
! IF (dd_max .ge. dd_his_max) then
! dd_his_max=dd_max
! ENDIF
! Output the average estrain in parts
ave_estran=pzero
! DO ip=1,NPart
! ave_estran(1) = ave_estran(1)+gestran(1,ip)*VolFrac(ip)
! ave_estran(2) = ave_estran(2)+gestran(2,ip)*VolFrac(ip)
! ave_estran(3) = ave_estran(3)+gestran(3,ip)*VolFrac(ip)
! ave_estran(4) = ave_estran(4)+gestran(4,ip)*VolFrac(ip)
! ave_estran(5) = ave_estran(5)+gestran(5,ip)*VolFrac(ip)
! ave_estran(6) = ave_estran(6)+gestran(6,ip)*VolFrac(ip)
! ENDDO
! Output the history max of the max spatial dd discrepancy
! IF (det_dd_max .ge. det_dd_max_his) then
! det_dd_max_his=det_dd_max
! ENDIF
!-->damin: comment this part out, incompatible with direct simulation
!***********************************************************************
call SaveStateVars(statev, nstatv, gstress, gestran, gkappa, &
ggamdot, grss, gmudot, gmu, geqvalues, gdd_for, gdd_revp, &
gdd_revn, gdd_stressf, gdd_rev0, gdd_deb, dd_all, dd_max, &
dd_his, dd_his_max, det_dd_max, ave_estran, det_dd_max_his, &
grainid)
! stresses and algorithmic moduli in abaqus format
CALL SaveStressModuli(stress, ddsdde, savg_ij, cavg_ijkl)
if(debug==1 .and. noel==1) write(File_e,*) 'kinc=',kinc
if(debug==1 .and. noel==1) write(File_e,'(6E20.10)') stress
if(debug==1 .and. noel==1) write(File_e,'(6E20.10)') strain
1000 FORMAT(/'*---time(1)---dtime---kinc---*'/, 7x, F16.10,F16.10,i8)
RETURN
END
!***********************************************************************
SUBROUTINE RecoverStateVars(statev, nstatv, gstress_n, gestran_n, &
gkappa_n, ggamdot, grss, gmudot,geqvalues, gmu ,gmu_n, &
gdd_for_n, gdd_for, gdd_revp_n, gdd_revp, gdd_revn_n, &
gdd_revn,gdd_stressf_n, gdd_stressf, gdd_rev0_n, gdd_rev0, &
gdd_deb_n, gdd_deb, grainid)
!***********************************************************************
USE NumType
USE FILEIO
Use Plapar
IMPLICIT NONE
INTEGER nstatv,kinc
REAL(KIND=8) :: statev(nstatv)
REAL(KIND=8) :: gstress_n (NVEC)
REAL(KIND=8) :: gestran_n (NVEC)
REAL(KIND=8) :: gkappa_n (MaxNumSlip)
REAL(KIND=8) :: gcrot_n (NSD, NSD)
REAL(KIND=8) :: grrot_n (NSD, NSD)
REAL(KIND=8) :: gcrot0 (NSD, NSD)
REAL(KIND=8) :: geqvalues (NEQVA)
REAL(KIND=8) :: gmu (NVEC)
REAL(KIND=8) :: gmu_n (NVEC)
REAL(KIND=8) :: ggamdot (MaxNumSlip)
REAL(KIND=8) :: grss (MaxNumSlip)
REAL(KIND=8) :: gmudot(NVEC)
REAL(KIND=8) :: gdd_for_n (MaxNumSlip)
REAL(KIND=8) :: gdd_revp_n (MaxNumSlip)
REAL(KIND=8) :: gdd_revn_n (MaxNumSlip)
REAL(KIND=8) :: gdd_stressf_n (MaxNumSlip)
REAL(KIND=8) :: gdd_rev0_n (MaxNumSlip)
REAL(KIND=8) :: gdd_deb_n (MaxNumSlip)
REAL(KIND=8) :: gdd_for (MaxNumSlip)
REAL(KIND=8) :: gdd_revp (MaxNumSlip)
REAL(KIND=8) :: gdd_revn (MaxNumSlip)
REAL(KIND=8) :: gdd_stressf (MaxNumSlip)
REAL(KIND=8) :: gdd_rev0 (MaxNumSlip)
REAL(KIND=8) :: gdd_deb (MaxNumSlip)
INTEGER varsPerPart, ip, dex, id, numslip, grainid
REAL(KIND=8) :: time(2)
!-----------------------------------------------------------------------
! recover state variables from abaqus vector array
dex=0
if (Plap%crystalID(grainid) .eq. kHCP) then
numslip=PhSlip(1)
else if (Plap%crystalID(grainid) .eq. kBCC) then
numslip=PhSlip(2)
else
call RunTimeError( FILE_O,'Error: crystalID exceed the existing lattice!')
end if
! stress
dex = dex + 1
gstress_n(1:NVec)=statev(dex:(dex+NVEC-1))
if (debug==1) then
write(FILE_E,*) 'gstress_n:'
write(FILE_E,'(6F12.5)') gstress_n(:)
endif
! estrain
dex = dex + NVEC
gestran_n(1:NVec)=statev(dex:(dex+NVEC-1))
if (debug==1) then
write(FILE_E,*) 'estran_n:'
write(FILE_E,'(6F12.5)') gestran_n(:)
end if
! kappa
dex = dex + NVEC
gkappa_n(1:numslip)=statev(dex:(dex+NumSlip-1))
if (debug==1) then
write(FILE_E,*) 'gkappa_n:'
write(FILE_E,'(6F12.5)') gkappa_n(:)
endif
! VonMises, Shearate, gamtot
dex = dex + numslip
statev(dex:(dex+NEQVA/2-1))=pzero
if (debug==1) then
write(FILE_E,*) '*----Initial VonMises, Shearate, gamtot in ', ip, '-th part----'
write(FILE_E,'(6F12.5)') statev(dex:(dex+NEQVA/2-1))
endif
! gammadot
dex = dex + NEQVA/2
ggamdot(1:numslip)=statev(dex:(dex+numslip-1))
if (debug==1) then
write(FILE_E,*) 'ggamdot:'
write(FILE_E,'(6F12.5)') ggamdot(:)
endif
! mu
dex = dex + NumSlip
gmu_n(1:NVec)=statev(dex:(dex+NVEC-1))
if (debug==1) then
write(FILE_E,*) 'gmu_n:'
write(FILE_E,'(6F12.5)') gmu_n(:)
endif
! rss
dex = dex + NVEC
grss(1:numslip)=statev(dex:(dex+numslip-1))
if (debug==1) then
write(FILE_E,*) 'grss:'
write(FILE_E,'(6F12.5)') grss(:)
endif
! mudot
dex = dex + NumSlip
gmudot(1:NVec)=statev(dex:(dex+NVEC-1))
if (debug==1) then
write(FILE_E,*) 'gmu_n:'
write(FILE_E,'(6F12.5)') gmudot(:)
endif
! dd_for
dex = dex + NVEC
gdd_for_n(1:numslip)=statev(dex:(dex+numslip-1))
if (debug==1) then
write(FILE_E,*) 'dd_for:'
write(FILE_E,*) gdd_for_n(:)
endif
! dd_revp
dex = dex + NumSlip
gdd_revp_n(1:numslip)=statev(dex:(dex+numslip-1))
if (debug==1) then
write(FILE_E,*) 'dd_revp:'
write(FILE_E,'(6F12.5)') gdd_revp_n(:)
endif
! dd_revn
dex = dex + NumSlip
gdd_revn_n(1:numslip)=statev(dex:(dex+numslip-1))
if (debug==1) then
write(FILE_E,*) 'dd_revn:'
write(FILE_E,'(6F12.5)') gdd_revn_n(:)
endif
! stressf
dex = dex + NumSlip
gdd_stressf_n(1:numslip)=statev(dex:(dex+numslip-1))
if (debug==1) then
write(FILE_E,*) 'stressf:'
write(FILE_E,'(6F12.5)') gdd_stressf_n(:)
endif
! dd_rev0
dex = dex + NumSlip
gdd_rev0_n(1:numslip)=statev(dex:(dex+numslip-1))
if (debug==1) then
write(FILE_E,*) 'dd_rev0:'
write(FILE_E,*) gdd_rev0_n(:)
endif
! dd_deb
dex = dex + NumSlip
gdd_deb_n(1:numslip)=statev(dex:(dex+numslip-1))
if (debug==1) then
write(FILE_E,*) 'dd_deb:'
write(FILE_E,*) gdd_deb_n(:)
endif
RETURN
END
!***********************************************************************
SUBROUTINE FetchCrystalVariablesAtPart(gstress_n, gestran_n, gkappa_n, &
geqvalues,stress_n, estran_n, kappa_n, eqvalues, mu_n, &
gmu_n, dd_for_n, gdd_for_n, dd_revp_n, gdd_revp_n, &
dd_revn_n, gdd_revn_n, dd_stressf_n, gdd_stressf_n, &
dd_rev0_n, gdd_rev0_n, dd_deb_n, gdd_deb_n)
!***********************************************************************
USE NumType
USE FILEIO
IMPLICIT NONE
REAL(KIND=8) :: gstress_n (NVEC)
REAL(KIND=8) :: gestran_n (NVEC)
REAL(KIND=8) :: gkappa_n (MaxNumSlip)
REAL(KIND=8) :: geqvalues (NEQVA)
REAL(KIND=8) :: gmu_n (NVEC)
REAL(KIND=8) :: gdd_for_n (MaxNumSlip)
REAL(KIND=8) :: gdd_revp_n (MaxNumSlip)
REAL(KIND=8) :: gdd_revn_n (MaxNumSlip)
REAL(KIND=8) :: gdd_stressf_n (MaxNumSlip)
REAL(KIND=8) :: gdd_rev0_n (MaxNumSlip)
REAL(KIND=8) :: gdd_deb_n (MaxNumSlip)
REAL(KIND=8) :: dd_for_n (MaxNumSlip)
REAL(KIND=8) :: dd_revp_n (MaxNumSlip)
REAL(KIND=8) :: dd_revn_n (MaxNumSlip)
REAL(KIND=8) :: dd_stressf_n (MaxNumSlip)
REAL(KIND=8) :: dd_rev0_n (MaxNumSlip)
REAL(KIND=8) :: dd_deb_n (MaxNumSlip)
REAL(KIND=8) :: stress_n(NVEC), estran_n(NVEC), kappa_n(MaxNumSlip),mu_n(NVEC)
REAL(KIND=8) :: eqvalues(NEQVA)
!-----------------------------------------------------------------------
stress_n= gstress_n
estran_n= gestran_n
kappa_n = gkappa_n
mu_n=gmu_n
dd_for_n=gdd_for_n
dd_revp_n=gdd_revp_n
dd_revn_n=gdd_revn_n
dd_stressf_n=gdd_stressf_n
dd_rev0_n=gdd_rev0_n
dd_deb_n=gdd_deb_n
eqvalues(kMISES_n) = geqvalues(kMISES_n)
eqvalues(kSHRATE_n) = geqvalues(kSHRATE_n)
eqvalues(kGAMTOT_n) = geqvalues(kGAMTOT_n)
RETURN
END
!***********************************************************************
SUBROUTINE ResetCrystalQnts(stress, estran, kappa, statev, eqvalues, &
gamdot, rss, mudot, mu, mu_n, s_ij, c_ijkl, stress_n, &
estran_n, kappa_n, dtime, VecM, VecS, ZVec, dd_for_n, &
dd_for, dd_revp_n, dd_revp, dd_revn_n, dd_revn, &
dd_stressf_n, dd_stressf, dd_rev0_n, dd_rev0, dd_deb_n, &
dd_deb, argmin, argmax, grainid )
!***********************************************************************
USE FILEIO
USE NumType
use PlaPar
IMPLICIT NONE
INTEGER :: is, grainid, numslip
REAL(KIND=8) :: dtime, argmin,argmax
REAL(KIND=8) :: stress(NVEC), estran(NVEC), kappa(MaxNumSlip)
REAL(KIND=8) :: statev(NSTAV)
REAL(KIND=8) :: s_ij(NSD,NSD), c_ijkl(NVEC,NVEC)
REAL(KIND=8) :: stress_n(NVEC), estran_n(NVEC), kappa_n(MaxNumSlip)
REAL(KIND=8) :: gamdot(MaxNumSlip), mu(NVEC),mu_n(NVEC)
REAL(KIND=8) :: rss(MaxNumSlip), crss(MaxNumSlip)
REAL(KIND=8) :: InnerProductVec, SSKineticEqn
REAL(KIND=8) :: mudot(NVEC)
REAL(KIND=8) :: eqvalues(NEQVA)
REAL(KIND=8) :: dd_for (MaxNumSlip)
REAL(KIND=8) :: dd_revp (MaxNumSlip)
REAL(KIND=8) :: dd_revn (MaxNumSlip)
REAL(KIND=8) :: dd_stressf (MaxNumSlip)
REAL(KIND=8) :: dd_rev0 (MaxNumSlip)
REAL(KIND=8) :: dd_deb (MaxNumSlip)
REAL(KIND=8) :: dd_for_n (MaxNumSlip)
REAL(KIND=8) :: dd_revp_n (MaxNumSlip)
REAL(KIND=8) :: dd_revn_n (MaxNumSlip)
REAL(KIND=8) :: dd_stressf_n (MaxNumSlip)
REAL(KIND=8) :: dd_rev0_n (MaxNumSlip)
REAL(KIND=8) :: dd_deb_n (MaxNumSlip)
REAL(KIND = 8) :: VecM(3,MaxNumSlip)
REAL(KIND = 8) :: VecS(3,MaxNumSlip)
REAL(KIND = 8) :: ZTen(3,3,MaxNumSlip),ZVec(6,MaxNumSlip)
!-----------------------------------------------------------------------
stress=stress_n
estran=estran_n
kappa = kappa_n
dd_for=dd_for_n
dd_revp=dd_revp_n
dd_revn=dd_revn_n
dd_stressf=dd_stressf_n
dd_rev0=dd_rev0_n
dd_deb=dd_deb_n
eqvalues(kMISES) = eqvalues(kMISES_n)
eqvalues(kSHRATE) = eqvalues(kSHRATE_n)
eqvalues(kGAMTOT) = eqvalues(kGAMTOT_n)
! Cauchy stress (tensor form)
CALL Vec6x1ToMat3x3Symm(stress(1), s_ij(1,1), NSD)
mudot(:)=pzero
if (Plap%crystalID(grainid) .eq. kHCP) then
numslip=PhSlip(1)
else if (Plap%crystalID(grainid) .eq. kBCC) then
numslip=PhSlip(2)
else
call RunTimeError( FILE_O,'Error: crystalID exceed the existing lattice!')
end if
DO is = 1, numslip
! Resolve shear stresses
rss(is) = InnerProductVec(stress(1), ZVec(1,is), NVEC)
crss(is)= kappa(is)
! Shear strain rate
gamdot(is) = SSKineticEqn(rss(is),crss(is),kDDGAM, argmin, argmax,is,grainid)
mudot(:)=mudot(:)+gamdot(is)*ZVec(:, is)
ENDDO
! consistent tangent
!CALL PlasticModuli(c_ijkl, stress, kappa, dtime, VecM, VecS, ZVec)
RETURN
END
!***********************************************************************
SUBROUTINE SaveCrystalVariablesAtPart(stress, estran, kappa, gamdot, &
rss, mudot, eqvalues, gstress, gestran, gkappa, gmu, mu, &
ggamdot, grss, gmudot, geqvalues, dd_for, gdd_for, dd_revp, &
gdd_revp, dd_revn, gdd_revn, dd_stressf, gdd_stressf, &
dd_rev0, gdd_rev0, dd_deb, gdd_deb)
!***********************************************************************
Use NumType
USE FILEIO
IMPLICIT NONE
INTEGER :: igrn, iqpt, ielem,ip
REAL(KIND=8) :: gstress (NVEC)
REAL(KIND=8) :: gestran (NVEC)
REAL(KIND=8) :: gkappa (MaxNumSlip)
REAL(KIND=8) :: eqvalues (NEQVA)
REAL(KIND=8) :: geqvalues(NEQVA)
REAL(KIND=8) :: gmu (NVEC)
REAL(KIND=8) :: ggamdot (MaxNumSlip)
REAL(KIND=8) :: grss (MaxNumSlip)
REAL(KIND=8) :: stress(NVEC), estran(NVEC), kappa(MaxNumSlip)
REAL(KIND=8) :: gamdot(MaxNumSlip),mu(NVEC), rss(MaxNumSlip)
REAL(KIND=8) :: gmudot(NVEC),mudot(NVEC)
REAL(KIND=8) :: gdd_for (MaxNumSlip)
REAL(KIND=8) :: gdd_revp (MaxNumSlip)
REAL(KIND=8) :: gdd_revn (MaxNumSlip)
REAL(KIND=8) :: gdd_stressf (MaxNumSlip)
REAL(KIND=8) :: gdd_rev0 (MaxNumSlip)
REAL(KIND=8) :: gdd_deb (MaxNumSlip)
REAL(KIND=8) :: dd_for (MaxNumSlip)
REAL(KIND=8) :: dd_revp (MaxNumSlip)
REAL(KIND=8) :: dd_revn (MaxNumSlip)
REAL(KIND=8) :: dd_stressf (MaxNumSlip)
REAL(KIND=8) :: dd_rev0 (MaxNumSlip)
REAL(KIND=8) :: dd_deb (MaxNumSlip)
!-----------------------------------------------------------------------
gstress=stress
gestran=estran
gkappa =kappa
ggamdot=gamdot
grss =rss
gmudot =mudot
gmu=mu
gdd_for= dd_for
gdd_revp= dd_revp
gdd_revn= dd_revn
gdd_stressf= dd_stressf
gdd_rev0= dd_rev0
gdd_deb= dd_deb
geqvalues(kMISES) = eqvalues(kMISES)
geqvalues(kSHRATE) = eqvalues(kSHRATE)
geqvalues(kGAMTOT) = eqvalues(kGAMTOT)
RETURN
END
!***********************************************************************
SUBROUTINE SaveStateVars(statev, nstatv, gstress, gestran, gkappa, &
ggamdot, grss, gmudot, gmu, geqvalues, gdd_for, gdd_revp, &
gdd_revn, gdd_stressf, gdd_rev0, gdd_deb, dd_all, dd_max, &
dd_his, dd_his_max, det_dd_max, ave_estran, det_dd_max_his, &
grainid )
!***********************************************************************
USE NumType
USE FILEIO
Use Plapar
IMPLICIT NONE
INTEGER :: nstatv
REAL(KIND=8) :: statev(nstatv)
REAL(KIND=8) :: gstress (NVEC)
REAL(KIND=8) :: gestran (NVEC)
REAL(KIND=8) :: gkappa (MaxNumSlip)
REAL(KIND=8) :: gmu (NVEC)
REAL(KIND=8) :: ggamdot (MaxNumSlip)
REAL(KIND=8) :: grss (MaxNumSlip)
REAL(KIND=8) :: gdd_for (MaxNumSlip)
REAL(KIND=8) :: gdd_revp (MaxNumSlip)
REAL(KIND=8) :: gdd_revn (MaxNumSlip)
REAL(KIND=8) :: gdd_stressf (MaxNumSlip)
REAL(KIND=8) :: gdd_rev0 (MaxNumSlip)
REAL(KIND=8) :: gdd_deb (MaxNumSlip)
REAL(KIND=8) :: dd_all
REAL(KIND=8) :: dd_max
REAL(KIND=8) :: dd_his
REAL(KIND=8) :: dd_his_max
REAL(KIND=8) :: det_dd_max
REAL(KIND=8) :: ave_estran (NVEC)
REAL(KIND=8) :: det_dd_max_his
REAL(KIND=8) :: geqvalues (NEQVA)
REAL(KIND=8) :: gmudot(NVEC)
INTEGER :: grainid, dex, id ,numslip
!-----------------------------------------------------------------------
! save state variables in abaqus vector array
ave_estran=0.d0
dex=0
if (Plap%crystalID(grainid) .eq. kHCP) then
numslip=PhSlip(1)
else if (Plap%crystalID(grainid) .eq. kBCC) then
numslip=PhSlip(2)
else
call RunTimeError(FILE_O, 'Error: crystalID exceed the existing lattice!')
end if
! stress
dex = dex + 1
statev(dex:(dex+NVEC-1))=gstress(1:NVEC)
if (debug==1) then
Write(FILE_E,*) 'initial estimate gstress='
Write(FILE_E,*) gstress
end if
! estrain
dex = dex + NVEC
statev(dex:(dex+NVEC-1))=gestran(1:NVEC)
! kappa
dex = dex + NVEC
statev(dex:(dex+NumSlip-1))=gkappa(1:numslip)
! VonMises, Shearate, gamtot
dex = dex + numslip
statev(dex:(dex+NEQVA/2-1))=geqvalues ((NEQVA/2+1):NEQVA)
! gamdot
dex = dex + NEQVA/2
statev(dex:(dex+numslip-1))=ggamdot(1:numslip)
! average estran: this value is wrong
dex = dex + NumSlip
statev(dex:(dex+NVEC-1))=gmu(1:NVEC)
ave_estran=ave_estran+gmu(1:NVEC)!*VolFrac(grainid)
! rss
dex = dex + NVEC
statev(dex:(dex+numslip-1))=grss(1:numslip)
! mudot
dex = dex + NumSlip
statev(dex:(dex+NVEC-1))=gmudot(1:NVEC)
! dd_for
dex = dex + NVEC
statev(dex:(dex+numslip-1))=gdd_for(1:numslip)
! dd_revp
dex = dex + NumSlip
statev(dex:(dex+numslip-1))=gdd_revp(1:numslip)
! dd_revn
dex = dex + NumSlip
statev(dex:(dex+numslip-1))=gdd_revn(1:numslip)
! stressf
dex = dex + NumSlip
statev(dex:(dex+numslip-1))=gdd_stressf(1:numslip)
! dd_rev0
dex = dex + NumSlip
statev(dex:(dex+numslip-1))=gdd_rev0(1:numslip)
! dd_deb
dex = dex + NumSlip
statev(dex:(dex+numslip-1))=gdd_deb(1:numslip)
dex = dex + NumSlip-1
! dd in each part
statev(nstatv-11)=dd_all
! dd_his max in each part
statev(nstatv-10)=dd_his
statev(nstatv-9)=dd_max ! dd_max among all parts
statev(nstatv-8)=dd_his_max ! dd_history_max among all parts
statev(nstatv-7)=det_dd_max ! dd_det_max_DD in each part
statev((nstatv-6):(nstatv-1))=ave_estran(1:6) ! average plastic strain in each part
statev(nstatv)=det_dd_max_his ! dd_history_det_max_DD in each part
RETURN
END
!***********************************************************************
SUBROUTINE SaveStressModuli(stress, ddsdde, savg_ij, cavg_ijkl)
!***********************************************************************
USE NumType
USE FILEIO
IMPLICIT NONE
REAL(KIND=8) :: stress(NVEC), ddsdde(NVEC, NVEC)
REAL(KIND=8) :: savg_ij(NSD, NSD), cavg_ijkl(NVEC, NVEC)
!-----------------------------------------------------------------------
! stresses
stress=pzero
stress(1) = savg_ij(1,1)
stress(2) = savg_ij(2,2)
stress(3) = savg_ij(3,3)
stress(4) = savg_ij(1,2)
stress(5) = savg_ij(1,3)
stress(6) = savg_ij(2,3)
! algorithmic moduli
ddsdde=pzero
ddsdde=ddsdde+ cavg_ijkl(:,:)
RETURN
END