Newer
Older
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
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
! |
! ||===\\ |
! || \\ |
! || || //==\\ || || //==|| ||/==\\ |
! || || || || || || || || || || |
! || // || || || || || || || |
! ||===// \\==// \\==\\ \\==\\ || |
! |
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
! |
! VRM Feb. 2007 |
! |
!------------------------------------------------------------------------------|
!------------------------------------------------------------------------------|
subroutine vrm(J2d,J3d,compressibility,viscosity,pressure,plasticity_parameters,plasticity_type,is_plastic,flag_vrm_pb,straintot)
use constants
use invariants
!------------------------------------------------------------------------------|
!(((((((((((((((( Purpose of the routine ))))))))))))))))))))))))))))))))))))))
!------------------------------------------------------------------------------|
!This routine computes the rescaled viscosity in the element in the case it is
!above the failure criterion
!------------------------------------------------------------------------------|
!(((((((((((((((( declaration of the subroutine arguments ))))))))))))))))))))
!------------------------------------------------------------------------------|
implicit none
double precision :: J2d,J3d,viscosity,pressure,compressibility,straintot
character (len=8) plasticity_type
double precision plasticity_parameters(9)
logical :: is_plastic,flag_vrm_pb
!------------------------------------------------------------------------------|
!(((((((((((((((( declaration of the subroutine internal variables )))))))))))))
!------------------------------------------------------------------------------|
double precision :: yield,fail,theta,c,e2d,zeta
double precision :: sin_theta,cos_theta,sin_phi,cos_phi
double precision :: strain_soft_in,strain_soft_out,strain_soft_phi,phi,fact
double precision :: strain_soft_yield
!-------------------------------------------------------------------------------
!-------------------------------------------------------------------------------
flag_vrm_pb = .false.
select case (trim(plasticity_type))
case ('vM')
yield = plasticity_parameters(1)
! strain softening
strain_soft_in=plasticity_parameters(3)
if (strain_soft_in.gt.0.d0) then
strain_soft_out=plasticity_parameters(4)
strain_soft_yield=plasticity_parameters(2)
if (straintot.gt.strain_soft_in) then
yield = (strain_soft_yield-yield)/(strain_soft_out-strain_soft_in)*(straintot-strain_soft_out) + strain_soft_yield
yield = max(yield, strain_soft_yield)
end if
end if
e2d = sqrt(J2d)
fail = 2.d0*viscosity*e2d-yield
if (fail.gt.0.d0) then
viscosity=0.5d0*yield/e2d
is_plastic=.true.
else
is_plastic=.false.
end if
case ('DPI','DPII')
!yield = abs(plasticity_parameters(9)-plasticity_parameters(8)*pressure)
! strain softening
phi=plasticity_parameters(1)
strain_soft_in=plasticity_parameters(3)
if (strain_soft_in.gt.0.d0) then
strain_soft_out=plasticity_parameters(4)
strain_soft_phi=plasticity_parameters(5)/180.d0*3.141592654d0
if (straintot.gt.strain_soft_in) then
fact=(straintot-strain_soft_in)/(strain_soft_out-strain_soft_in)
fact=min(fact,1.d0)
phi=phi+(strain_soft_phi-phi)*fact
endif
endif
yield = max(plasticity_parameters(2)/compressibility,plasticity_parameters(2)+tan(phi)*pressure)
e2d = sqrt(J2d)
fail = 2.d0*viscosity*e2d-yield
if (fail.gt.0.d0) then
viscosity=0.5d0*yield/e2d
is_plastic=.true.
else
is_plastic=.false.
end if
case ('DPIII')
e2d = sqrt(J2d)
!yield = plasticity_parameters(2)+plasticity_parameters(1)*pressure
yield = max(plasticity_parameters(2)/compressibility,plasticity_parameters(2)+plasticity_parameters(1)*pressure)
fail = 2.d0*viscosity*e2d-yield
if (fail.gt.0.d0) then
viscosity=0.5d0*yield/e2d
is_plastic=.true.
else
is_plastic=.false.
end if
case ('MC')
e2d = sqrt(J2d)
theta = lode_angle(J2d,J3d)
sin_theta = sin(theta)
cos_theta = cos(theta)
sin_phi = plasticity_parameters(8)
cos_phi = plasticity_parameters(9)
c = plasticity_parameters(2)
! strain softening
strain_soft_in=plasticity_parameters(3)
if (strain_soft_in.gt.0.d0) then
strain_soft_out=plasticity_parameters(4)
strain_soft_phi=plasticity_parameters(5)/180.d0*3.141592654d0
if (straintot.gt.strain_soft_in) then
phi=atan2(plasticity_parameters(8),plasticity_parameters(9))
fact=(straintot-strain_soft_in)/(strain_soft_out-strain_soft_in)
fact=min(fact,1.d0)
phi=phi+(strain_soft_phi-phi)*fact
sin_phi=sin(phi)
cos_phi=cos(phi)
endif
endif
zeta = cos_theta-sin_phi*sin_theta/sqrt3
! modified by Jean
! fail = -pressure*sin_phi+e2d*zeta-c*cos_phi
fail = pressure*sin_phi+2.d0*viscosity*e2d*zeta-c*cos_phi
if (fail.gt.0.d0) then
viscosity=0.5d0/e2d*(c*cos_phi-pressure*sin_phi)/zeta
!if (viscosity<0.d0) then
! viscosity=0.d0
! !print *,viscosity,c*cos_phi,pressure*sin_phi
!end if
is_plastic=.true.
else
is_plastic=.false.
end if
case ('Tresca')
e2d = sqrt(J2d)
theta = lode_angle(J2d,J3d)
yield = plasticity_parameters(1)
cos_theta = cos(theta)
fail = 2.d0*viscosity*e2d*cos_theta-yield
if (fail.gt.0.d0) then
viscosity=0.5d0*yield/e2d/cos_theta
is_plastic=.true.
else
is_plastic=.false.
end if
case default
print *,plasticity_type
call stop_run('error in vrm$')
flag_vrm_pb = .true.
end select
return
end
!-------------------------------------------------------------------------------
!-------------------------------------------------------------------------------
subroutine compute_plastic_params (params,mat)
use definitions
use constants
implicit none
type (parameters) params
type (material) :: mat(0:params%nmat)
double precision c,k,alpha,phi
integer i,iproc,nproc,ierr
INCLUDE 'mpif.h'
call mpi_comm_size (mpi_comm_world,nproc,ierr)
call mpi_comm_rank (mpi_comm_world,iproc,ierr)
do i=0,params%nmat
select case (trim(mat(i)%plasticity_type))
case ('No')
case ('vM')
case ('DPI')
phi = mat(i)%plasticity_parameters(1)*pi/180.d0
c = mat(i)%plasticity_parameters(2)
alpha = 2.d0*sin(phi) /sqrt(3.d0)/(3.d0-sin(phi))
k = 6.d0*c*cos(phi)/sqrt(3.d0)/(3.d0-sin(phi))
mat(i)%plasticity_parameters(8)=alpha
mat(i)%plasticity_parameters(9)=k
if (iproc.eq.0) then
write(8,*) 'material ',i
write(8,*) ' phi =',phi
write(8,*) ' c =',c
write(8,*) ' alpha=',alpha
write(8,*) ' k =',k
end if
case ('DPII')
phi = mat(i)%plasticity_parameters(1)*pi/180.d0
c = mat(i)%plasticity_parameters(2)
alpha = tan(phi)/sqrt(9.d0+12.d0*(tan(phi))**2)
k = 3.d0*c /sqrt(9.d0+12.d0*(tan(phi))**2)
mat(i)%plasticity_parameters(8)=alpha
mat(i)%plasticity_parameters(9)=k
case ('MC')
phi = mat(i)%plasticity_parameters(1)*pi/180.d0
mat(i)%plasticity_parameters(8)=sin(phi)
mat(i)%plasticity_parameters(9)=cos(phi)
case ('Tresca')
case ('DPIII')
case default
call stop_run('error in compute_plastic_params$')
end select
if (iproc.eq.0) then
write(8,*) 'material ',i,'plasticity_type ',mat(i)%plasticity_type
write(8,*) 'plasticity_parameters',mat(i)%plasticity_parameters(:)
end if
end do
return
end