-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathmodel.py
368 lines (305 loc) · 15.5 KB
/
model.py
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
import numpy as np
from scipy import sparse
import random
from Stokeselvis import return_sparse_matrix_Stokes
from interpolate import interpolate, interpolate2m, interpolate_harmonic, fill_nans
from largetime import Time
import matplotlib.pyplot as plt
average = lambda x: (x[:-1,:-1] + x[1:,:-1] +x[1:,1:] +x[:-1,1:])/4.0
class Model(object):
def __init__(self, parameters):
defaults = {
'T' : Time(0),
'step' : 0,
'width' : None,
'height' : None,
'j_res' : None,
'i_res' : None,
'gx_0' : None,
'gy_0' : None,
'p0cell' : 0,
'right_bound' : None,
'left_bound' : None,
'top_bound' : None,
'bottom_bound' : None,
'stress_changes' : 'simple', # '2nd order', 'subgird'
'subgrid_relaxation' : .5,
'advect_scheme' : 'simple', # 'Runge-Kutta 2nd order'
'velocity_multiplier' : 10**-7,
}
self.markers = {'mxx' : None,
'myy' : None,
'm_cat' : None,
'm_rho' : None,
'm_eta' : None,
'm_C' : None,
'm_mu' : None,
'm_sinphi' : None,
'm_s_xx' : None,
'm_s_xy' : None,
'm_e_xx' : None,
'm_e_xy' : None,
'm_P' : None,
'markers_index_list' : None,
'moving_points_index_list' : None,
}
self.__dict__.update(defaults)
self.__dict__.update(self.markers)
for key in self.__dict__:
try:
self.__dict__[key] = parameters[key]
except KeyError:
pass
if self.__dict__[key] is None:
raise ValueError (f'{key} parameter must be set')
self.width = int(self.width)
self.height = int(self.height)
self.i_res = int(self.i_res)
self.j_res = int(self.j_res)
self.dx = self.width / (self.j_res-1)
self.dy = self.height /( self.i_res-1)
self.dx2, self.dy2 = self.dx**2, self.dy**2
self.kbond = 4*self.m_eta.min()/((self.dx+self.dy)**2)
self.kcont = 2*self.m_eta.min()/(self.dx+self.dy)
def make_step(self, maxT = 1):
dt_min=1e+10
mxx = self.mxx
myy = self.myy
m_cat = self.m_cat
m_rho = self.m_rho
m_eta = self.m_eta
m_mu = self.m_mu
m_C = self.m_C
m_sinphi = self.m_sinphi
m_s_xx = self.m_s_xx
m_s_xy = self.m_s_xy
m_e_xx = self.m_e_xx
m_e_xy = self.m_e_xy
m_P = self.m_P
m_s_ii_old = np.zeros(np.shape(mxx))
m_e_ii_old = np.zeros(np.shape(mxx))
i_res, j_res = self.i_res, self.j_res
dx, dy = self.dx, self.dy
gx_0, gy_0 = self.gx_0, self.gy_0
kbond, kcont = self.kbond, self.kcont
p0cell = self.p0cell
eta_min = 1e+10
eta_max = 1e+34
T = self.T
step = self.step
if not T : T = Time(0)
if not step: step = -1
dt = 0
plastic_iterations = 1
while T < maxT:
dt = max(dt, dt_min)
for i in range(plastic_iterations):
plastic_yield = False
m_xelvis = m_eta/(m_mu*dt + m_eta)
m_s_xx_new = m_s_xx*m_xelvis + 2*m_eta*m_e_xx*(1-m_xelvis)
m_s_xy_new = m_s_xy*m_xelvis + 2*m_eta*m_e_xy*(1-m_xelvis)
m_s_ii_new = (m_s_xx_new**2 + m_s_xy_new**2)**.5
m_s_ii_yield = m_C + m_sinphi*m_P
mask = m_s_ii_yield <0
m_s_ii_yield[mask] = 0
ym = m_s_ii_yield > m_s_ii_new # yielding mask
check_for_plastic = False
#print("step")
if check_for_plastic and np.all(~ym):
#ym = m_s_ii_yield < m_s_ii_new # yielding mask
#if not np.all(ym):
plastic_yield = True
print("plastic yield!")
m_e_ii_old[ym] = np.sqrt(m_e_xx[ym]**2 + m_e_xy[ym]**2)
m_eta[ym] = m_s_ii_yield[ym] /2/m_e_ii_old[ym]
m_s_ii_old[ym] = np.sqrt(m_s_xx[ym]**2 + m_s_xy[ym]**2)
m_s_xx[ym] = m_s_xx[ym]*m_s_ii_yield[ym]/m_s_ii_old[ym]
m_s_xy[ym] = m_s_xy[ym]*m_s_ii_yield[ym]/m_s_ii_old[ym]
# eta_min < eta < eta_max
# m_eta[m_eta<eta_min] = eta_min
# m_eta[m_eta>eta_max] = eta_max
eta_s, rho, so_xy = interpolate(mxx,myy,i_res,j_res, (m_eta, m_rho, m_s_xy))
eta_n, so_xx = interpolate(mxx+.5,myy+.5,i_res,j_res, (m_eta, m_s_xx))
mu_s = interpolate_harmonic(mxx,myy,i_res,j_res, m_mu )
mu_n = interpolate_harmonic(mxx+.5,myy+.5,i_res,j_res, m_mu )
eta_n = interpolate_harmonic(mxx+.5,myy+.5,i_res,j_res, m_eta )
eta_s = interpolate_harmonic(mxx,myy,i_res,j_res, m_eta )
# rho = interpolate_harmonic(mxx,myy,i_res,j_res, m_rho )
# so_xy = interpolate_harmonic(mxx,myy,i_res,j_res, m_s_xx )
# so_xx = interpolate_harmonic(mxx+.5,myy+.5,i_res,j_res, m_s_xx )
#Check if we have nans
if np.isnan(eta_s).any(): fill_nans(eta_s)
if np.isnan(eta_n).any(): fill_nans(eta_n)
if np.isnan(rho).any(): fill_nans(rho)
if np.isnan(mu_s).any(): fill_nans(mu_s)
if np.isnan(mu_n).any(): fill_nans(mu_n)
if np.isnan(so_xx).any(): fill_nans(so_xx)
if np.isnan(so_xy).any(): fill_nans(so_xy)
# compute viscoelastic (numerical) viscosity and stress
eta_s0 = eta_s
eta_n0 = eta_n
so_xy0 = so_xy
so_xx0 = so_xx
xelvis_s = eta_s/(eta_s+dt*mu_s)
xelvis_n = eta_n/(eta_n+dt*mu_n)
eta_s = eta_s*(1-xelvis_s)
eta_n = eta_n*(1-xelvis_n)
so_xy = so_xy*xelvis_s
so_xx = so_xx*xelvis_n
gx_0, gy_0 = self.gx_0, self.gy_0
Vbound = {}
for index, VxVy in self.moving_points_index_list:
x,y = mxx[index][0], myy[index][0]
# i,j = int(round(y)), int(round(x))
i,j = int(y), int(x)
# velocities defined as a function
if callable(VxVy):
Vx, Vy = VxVy(step=step, time=T)
else:
Vx, Vy = VxVy
Vx *= self.velocity_multiplier
Vy *= self.velocity_multiplier
Vbound[(i,j)] = [Vx,Vy]
print(Vbound)
eta_s_ = np.copy(eta_s)
eta_n_ = np.copy(eta_n)
rho_ = np.copy(rho)
so_xx_ = np.copy(so_xx)
so_xy_ = np.copy(so_xy)
Stokes_sparse, vector = return_sparse_matrix_Stokes(j_res, i_res, dx, dy,
eta_s, eta_n, rho, gx_0, gy_0, so_xx, so_xy, kbond, kcont, p0cell,
lower_boundary=self.bottom_bound, upper_boundary=self.top_bound,
right_boundary=self.right_bound, left_boundary=self.left_bound,
Vbound=Vbound)
Stokes_solve = sparse.linalg.spsolve(Stokes_sparse, vector)
P = Stokes_solve[::3].reshape((i_res),(j_res))
Vx = Stokes_solve[1::3].reshape((i_res),(j_res))
Vy = Stokes_solve[2::3].reshape((i_res),(j_res))
P *= kcont
Vx_max = np.abs(Vx).max()
Vy_max = np.abs(Vy).max()
dtx = 0.2 * dx / Vx_max
dty = 0.2 * dy / Vy_max
dt = min(dtx,dty)
#dVx, dVy = Vx[1:, :-1] - Vx[:-1, :-1], Vy[ :-1, 1:] - Vy[:-1, :-1]
#dVx_, dVy_ = Vx[ :, 1:-1] - Vx[:, :-2], Vy[1:-1, :] - Vy[:-2, :]
dVx_dx = (Vx[:-1,1:] - Vx[:-1,:-1])/dx
dVy_dy = (Vy[1:,:-1] - Vy[:-1,:-1])/dy
#dVx_dy = (Vx[:-1,1:] - Vx[:-1,:-1])/dy
#dVy_dx = (Vy[1:,:-1] - Vy[:-1,:-1])/dx
dVx_dy = (Vx[1:-1,:] - Vx[:-2,:])/dy
dVy_dx = (Vy[:,1:-1] - Vy[:,:-2])/dx
#concatenating array to add additional rows and columns. Just copy remaning columns
# abc ad aadd
# abc -> abc be -> bbee
# def def cf ccff
# def
dVx_dy = np.append(dVx_dy[0,:][None,:],dVx_dy,axis=0)
dVx_dy = np.append(dVx_dy,dVx_dy[-1,:][None,:],axis=0)
dVy_dx = np.append(dVy_dx[:,0][:,None],dVy_dx,axis=1)
dVy_dx = np.append(dVy_dx,dVy_dx[:,-1][:,None],axis=1)
e_xx = dVx_dx # strain rate
e_xy = .5 * (dVx_dy + dVy_dx)
s_xx = (1-xelvis_n[1:,1:])*2*eta_n0[1:,1:]*e_xx + xelvis_n[1:,1:]*so_xx0[1:,1:]
s_xy = (1-xelvis_s)*2*eta_s0*e_xy + xelvis_s*so_xy0
s_ii = (s_xx**2 + average(s_xy**2))**.5
e_ii = (e_xx**2 + average(e_xy**2))**.5
# ds_xx = s_xx - average(so_xx0)
ds_xx = s_xx - so_xx0[1:,1:]
ds_xy = s_xy - so_xy0
if not plastic_yield: break
m_Vx = interpolate2m(mxx , myy-.5, Vx[:-1,:])
m_Vy = interpolate2m(mxx-.5, myy , Vy[:,:-1])
m_P = interpolate2m(mxx-.5, myy-.5, P[1:,1:]) # tecnichaly, there must be +.5,+.5
# but since we slice P, indexing goes one item lower
w = 0.5*(dVy_dx - dVx_dy)
m_s_xx, m_s_xy = self.interpolate_stress_changes(mxx, myy,
m_s_xx, m_s_xy,
s_xx, s_xy,
ds_xx, ds_xy,
dt)
m_e_xx = interpolate2m(mxx-.5,myy-.5,e_xx)
m_e_xy = interpolate2m(mxx,myy,e_xy)
m_w = interpolate2m(mxx-.5 ,myy-.5 , w)
m_a = m_w * dt
# m_s_xx_ = m_s_xx - m_s_xy * 2 * m_a
# m_s_xy_ = m_s_xy + m_s_xx * 2 * m_a
# m_s_xx, m_s_xy = m_s_xx_, m_s_xy_
m_s_xx_ = m_s_xx*(np.cos(m_a)**2 - np.sin(m_a)**2 ) - m_s_xy*np.sin(2*m_a)
m_s_xy_ = m_s_xx*np.sin(2*m_a) + m_s_xy*np.cos(2*m_a)
self.advect(mxx, myy, m_Vx, m_Vy, Vx, Vy, dt)
T += dt
step +=1
parameters = {'T' : T,
'step' : step,
'eta_s' : eta_s,
'eta_n' : eta_n,
'mxx' : mxx,
'myy' : myy,
'm_cat' : m_cat,
'sii' : s_ii,
'eii' : e_ii,
'P' : P,
'Vx' : Vx,
'Vy' : Vy,
'e_xx' : e_xx,
'e_xy' : e_xy,
's_xx' : s_xx,
's_xy' : s_xy,
'xelvis_s' : xelvis_s,
'mu_n' : mu_n,
'mu_s' : mu_s,
'w' : w,
'markers_index_list' : self.markers_index_list,
'moving_points_index_list' : self.moving_points_index_list,
'eta_s_' : eta_s_,
'eta_n_' : eta_n_,
'rho_' : rho_,
'so_xx_' : so_xx_,
'so_xy_' : so_xy_,
}
yield parameters
def interpolate_stress_changes(self, mxx, myy, m_s_xx, m_s_xy, s_xx=None, s_xy=None, ds_xx=None, ds_xy=None, dt=None ):
if self.stress_changes == 'simple':
m_s_xx = interpolate2m(mxx-.5,myy-.5,s_xx)
m_s_xy = interpolate2m(mxx,myy,s_xy)
elif self.stress_changes == '2nd order':
print('2nd order')
m_ds_xx = interpolate2m(mxx-.5,myy-.5,ds_xx)
m_ds_xy = interpolate2m(mxx,myy,ds_xy)
m_s_xx = m_s_xx + m_ds_xx
m_s_xy = m_s_xy + m_ds_xy
elif self.stress_changes == 'subgrid':
print('subgrid')
d_ve = self.subgrid_relaxation
m_eta = self.m_eta
m_mu = self.m_mu
i_res, j_res = self.i_res, self.j_res
m_dt_maxwell = m_eta / m_mu
d_ = (-d_ve*dt/m_dt_maxwell).astype(np.float)
m_s_xx_nodes = interpolate2m(mxx-.5,myy-.5,s_xx)
m_s_xy_nodes = interpolate2m(mxx,myy,s_xy)
m_ds_xx_subgrid = (m_s_xx_nodes - m_s_xx)*(1-np.exp(d_))
m_ds_xy_subgrid = (m_s_xy_nodes - m_s_xy)*(1-np.exp(d_))
(ds_xx_subgrid,) = interpolate(mxx+.5,myy+.5,i_res,j_res, (m_ds_xx_subgrid, ))
(ds_xy_subgrid,) = interpolate(mxx,myy,i_res,j_res, (m_ds_xy_subgrid, ))
ds_xx_remaining = ds_xx - ds_xx_subgrid[1:,1:]
ds_xy_remaining = ds_xy - ds_xy_subgrid
m_ds_xx_remaining = interpolate2m(mxx-.5,myy-.5, ds_xx_remaining)
m_ds_xy_remaining = interpolate2m(mxx,myy, ds_xy_remaining)
m_s_xx += m_ds_xx_subgrid + m_ds_xx_remaining
m_s_xy += m_ds_xy_subgrid + m_ds_xy_remaining
return m_s_xx, m_s_xy
def advect(self, mxx, myy, m_Vx, m_Vy, Vx, Vy, dt ):
dx, dy = self.dx, self.dy
if self.advect_scheme=='simple':
mxx += m_Vx*dt/dx
myy += m_Vy*dt/dy
elif self.advect_scheme=='Runge-Kutta 2nd order':
mxx_b = mxx + 0.5*m_Vx*dt/dx
myy_b = myy + 0.5*m_Vy*dt/dy
m_Vx_eff = interpolate2m(mxx_b , myy_b-.5, Vx[:-1,:])
m_Vy_eff = interpolate2m(mxx_b-.5, myy_b , Vy[:,:-1])
mxx += m_Vx_eff*dt/dx
myy += m_Vy_eff*dt/dy
return mxx, myy