-
Notifications
You must be signed in to change notification settings - Fork 3
/
Copy pathStokes2D_vep_reg_IU.jl
167 lines (167 loc) · 7.14 KB
/
Stokes2D_vep_reg_IU.jl
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
using Plots,LinearAlgebra,Printf
# helper functions
@views av(A) = 0.25*(A[1:end-1,1:end-1].+A[2:end,1:end-1].+A[1:end-1,2:end].+A[2:end,2:end])
@views av_xa(A) = 0.5*(A[1:end-1,:].+A[2:end,:])
@views av_ya(A) = 0.5*(A[:,1:end-1].+A[:,2:end])
@views maxloc(A) = max.(A[1:end-2,1:end-2],A[1:end-2,2:end-1],A[1:end-2,3:end],
A[2:end-1,1:end-2],A[2:end-1,2:end-1],A[2:end-1,3:end],
A[3:end ,1:end-2],A[3:end ,2:end-1],A[3:end ,3:end])
@views bc2!(A) = begin A[1,:] = A[2,:]; A[end,:] = A[end-1,:]; A[:,1] = A[:,2]; A[:,end] = A[:,end-1] end
# main function
@views function Stokes2D_vep()
# phyiscs
lx,ly = 1.0,1.0
radi = 0.1
η0 = 1.0
η_reg = 1.2e-2
G0 = 1.0
Gi = 0.5G0
τ_y = 1.6
sinϕ = sind(30)
ebg = 1.0
dt = η0/G0/4.0
# numerics
nx,ny = 63,63
nt = 10
εnl = 1e-6
maxiter = 150max(nx,ny)
nchk = 2max(nx,ny)
Re = 5π
r = 1.0
CFL = 0.99/sqrt(2)
# preprocessing
dx,dy = lx/nx,ly/ny
max_lxy = max(lx,ly)
vpdτ = CFL*min(dx,dy)
xc,yc = LinRange(-(lx-dx)/2,(lx-dx)/2,nx),LinRange(-(ly-dy)/2,(ly-dy)/2,ny)
xv,yv = LinRange(-lx/2,lx/2,nx+1),LinRange(-ly/2,ly/2,ny+1)
# allocate arrays
Pr = zeros(nx ,ny )
τxx = zeros(nx ,ny )
τyy = zeros(nx ,ny )
τxy = zeros(nx+1,ny+1)
τxyc = zeros(nx ,ny )
τii = zeros(nx ,ny )
Eii = zeros(nx ,ny )
λ = zeros(nx ,ny )
F = zeros(nx ,ny )
Fchk = zeros(nx ,ny )
dQdTxx = zeros(nx ,ny )
dQdTyy = zeros(nx ,ny )
dQdTxy = zeros(nx ,ny )
τxx_o = zeros(nx ,ny )
τyy_o = zeros(nx ,ny )
τxy_o = zeros(nx+1,ny+1)
Vx = zeros(nx+1,ny )
Vy = zeros(nx ,ny+1)
dVx = zeros(nx-1,ny )
dVy = zeros(nx ,ny-1)
Rx = zeros(nx-1,ny )
Ry = zeros(nx ,ny-1)
∇V = zeros(nx ,ny )
ρg = zeros(nx ,ny )
Exx = zeros(nx ,ny )
Eyy = zeros(nx ,ny )
Exy = zeros(nx+1,ny+1)
Exx_e = zeros(nx ,ny )
Eyy_e = zeros(nx ,ny )
Exy_e = zeros(nx+1,ny+1)
Exx_τ = zeros(nx ,ny )
Eyy_τ = zeros(nx ,ny )
Exy_τ = zeros(nx+1,ny+1)
Exyc_τ = zeros(nx ,ny )
η_ve_τ = zeros(nx ,ny )
η_ve_τv = zeros(nx+1,ny+1)
η_ve = zeros(nx ,ny )
η_vem = zeros(nx ,ny )
η_vev = zeros(nx+1,ny+1)
η_vevm = zeros(nx+1,ny+1)
dτ_ρ = zeros(nx ,ny )
dτ_ρv = zeros(nx+1,ny+1)
Gdτ = zeros(nx ,ny )
Gdτv = zeros(nx+1,ny+1)
# init
Vx = [ ebg*x for x ∈ xv, _ ∈ yc ]
Vy = [-ebg*y for _ ∈ xc, y ∈ yv ]
η = fill(η0,nx,ny); ηv = fill(η0,nx+1,ny+1)
G = fill(G0,nx,ny); Gv = fill(G0,nx+1,ny+1)
@. G[xc^2 + yc'^2 < radi^2] = Gi
@. Gv[xv^2 + yv'^2 < radi^2] = Gi
η_e = G.*dt; η_ev = Gv.*dt
@. η_ve = 1.0/(1.0/η + 1.0/η_e)
@. η_vev = 1.0/(1.0/ηv + 1.0/η_ev)
# action
t = 0.0; evo_t = Float64[]; evo_τxx = Float64[]
for it = 1:nt
τxx_o .= τxx; τyy_o .= τyy; τxy_o .= τxy
err = 2εnl; iter = 0
while err > εnl && iter < maxiter
η_vem[2:end-1,2:end-1] .= maxloc(η_ve) ; bc2!(η_vem)
η_vevm[2:end-1,2:end-1] .= maxloc(η_vev); bc2!(η_vevm)
@. dτ_ρ = vpdτ*max_lxy/Re/η_vem
@. dτ_ρv = vpdτ*max_lxy/Re/η_vevm
@. Gdτ = vpdτ^2/dτ_ρ/(r+2.0)
@. Gdτv = vpdτ^2/dτ_ρv/(r+2.0)
@. η_ve_τ = 1.0/(1.0/η + 1.0/η_e + 1.0/Gdτ)
@. η_ve_τv = 1.0/(1.0/ηv + 1.0/η_ev + 1.0/Gdτv)
# pressure
∇V .= diff(Vx, dims=1)./dx .+ diff(Vy, dims=2)./dy
@. Pr -= r*Gdτ*∇V
# strain rates
Exx .= diff(Vx, dims=1)./dx
Eyy .= diff(Vy, dims=2)./dy
Exy[2:end-1,2:end-1] .= 0.5.*(diff(Vx[2:end-1,:], dims=2)./dy .+ diff(Vy[:,2:end-1], dims=1)./dx)
# viscoelastic strain rates
@. Exx_e = Exx + τxx_o/2.0/η_e
@. Eyy_e = Eyy + τyy_o/2.0/η_e
@. Exy_e = Exy + τxy_o/2.0/η_ev
# viscoelastic pseudo-transient strain rates
@. Exx_τ = Exx_e + τxx/2.0/Gdτ
@. Eyy_τ = Eyy_e + τyy/2.0/Gdτ
@. Exy_τ = Exy_e + τxy/2.0/Gdτv
# stress update
@. τxx = 2.0*η_ve_τ *Exx_τ
@. τyy = 2.0*η_ve_τ *Eyy_τ
@. τxy = 2.0*η_ve_τv*Exy_τ
@. τxyc = 2.0*η_ve_τ*Exyc_τ
# stress and strain rate invariants
@. τii = sqrt(0.5*(τxx^2 + τyy^2) + τxyc*τxyc)
@. Eii = sqrt(0.5*(Exx_τ^2 + Eyy_τ^2) + Exyc_τ*Exyc_τ)
# yield function
@. F = τii - τ_y - Pr.*sinϕ
@. λ = max(F,0.0)/(η_ve_τ + η_reg)
@. dQdTxx = 0.5*τxx /τii
@. dQdTyy = 0.5*τyy /τii
@. dQdTxy = τxyc/τii
# plastic correction
@. τxx = 2.0*η_ve_τ *(Exx_τ - λ*dQdTxx)
@. τyy = 2.0*η_ve_τ *(Eyy_τ - λ*dQdTyy)
@. τxyc = 2.0*η_ve_τ *(Exyc_τ - 0.5*λ*dQdTxy)
τxy[2:end-1,2:end-1] .= 2.0 .* η_ve_τv[2:end-1,2:end-1].*(Exy_τ[2:end-1,2:end-1] .- 0.5 .* av(λ.*dQdTxy))
@. τii = sqrt(0.5*(τxx^2 + τyy^2) + τxyc*τxyc)
@. Fchk = τii - τ_y - Pr*sinϕ - λ*η_reg
# velocity update
dVx .= av_xa(dτ_ρ) .* (.-diff(Pr, dims=1)./dx .+ diff(τxx, dims=1)./dx .+ diff(τxy[2:end-1,:], dims=2)./dy)
dVy .= av_ya(dτ_ρ) .* (.-diff(Pr, dims=2)./dy .+ diff(τyy, dims=2)./dy .+ diff(τxy[:,2:end-1], dims=1)./dx .+ av_ya(ρg))
@. Vx[2:end-1,:] = Vx[2:end-1,:] + dVx
@. Vy[:,2:end-1] = Vy[:,2:end-1] + dVy
if iter % nchk == 0
Rx .= .-diff(Pr, dims=1)./dx .+ diff(τxx, dims=1)./dx .+ diff(τxy[2:end-1,:], dims=2)./dy
Ry .= .-diff(Pr, dims=2)./dy .+ diff(τyy, dims=2)./dy .+ diff(τxy[:,2:end-1], dims=1)./dx .+ av_ya(ρg)
norm_Rx = norm(Rx)/sqrt(length(Rx)); norm_Ry = norm(Ry)/sqrt(length(Ry)); norm_∇V = norm(∇V)/sqrt(length(∇V))
err = maximum([norm_Rx, norm_Ry, norm_∇V])
@printf("it = %d, iter = %d, err = %1.2e norm[Rx=%1.2e, Ry=%1.2e, ∇V=%1.2e] (Fchk=%1.2e) \n", it, iter, err, norm_Rx, norm_Ry, norm_∇V, maximum(Fchk))
end
iter += 1
end
t += dt; push!(evo_t, t); push!(evo_τxx, maximum(τxx))
p1 = heatmap(xc,yc,τii',aspect_ratio=1,xlims=(-lx/2,lx/2),ylims=(-ly/2,ly/2),title="τii")
p2 = plot(evo_t, evo_τxx , legend=false, xlabel="time", ylabel="max(τxx)", linewidth=0, markershape=:circle, framestyle=:box, markersize=3)
plot!(evo_t, 2.0.*ebg.*η0.*(1.0.-exp.(.-evo_t.*G0./η0)), linewidth=2.0) # analytical solution for VE loading
plot!(evo_t, 2.0.*ebg.*η0.*ones(size(evo_t)), linewidth=2.0) # viscous flow stress
display(plot(p1,p2))
end
return
end
# action
Stokes2D_vep()