/
Navier_Stokes_inverse.py
184 lines (169 loc) · 6.73 KB
/
Navier_Stokes_inverse.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
"""Backend supported: tensorflow.compat.v1, tensorflow, pytorch, paddle
An inverse problem of the Navier-Stokes equation of incompressible flow around cylinder with Re=100
References: https://doi.org/10.1016/j.jcp.2018.10.045 Section 4.1.1
"""
import deepxde as dde
import numpy as np
import matplotlib.pyplot as plt
from scipy.io import loadmat
import re
# true values
C1true = 1.0
C2true = 0.01
# Load training data
def load_training_data(num):
data = loadmat("../dataset/cylinder_nektar_wake.mat")
U_star = data["U_star"] # N x 2 x T
P_star = data["p_star"] # N x T
t_star = data["t"] # T x 1
X_star = data["X_star"] # N x 2
N = X_star.shape[0]
T = t_star.shape[0]
# Rearrange Data
XX = np.tile(X_star[:, 0:1], (1, T)) # N x T
YY = np.tile(X_star[:, 1:2], (1, T)) # N x T
TT = np.tile(t_star, (1, N)).T # N x T
UU = U_star[:, 0, :] # N x T
VV = U_star[:, 1, :] # N x T
PP = P_star # N x T
x = XX.flatten()[:, None] # NT x 1
y = YY.flatten()[:, None] # NT x 1
t = TT.flatten()[:, None] # NT x 1
u = UU.flatten()[:, None] # NT x 1
v = VV.flatten()[:, None] # NT x 1
p = PP.flatten()[:, None] # NT x 1
# training domain: X × Y = [1, 8] × [−2, 2] and T = [0, 7]
data1 = np.concatenate([x, y, t, u, v, p], 1)
data2 = data1[:, :][data1[:, 2] <= 7]
data3 = data2[:, :][data2[:, 0] >= 1]
data4 = data3[:, :][data3[:, 0] <= 8]
data5 = data4[:, :][data4[:, 1] >= -2]
data_domain = data5[:, :][data5[:, 1] <= 2]
# choose number of training points: num =7000
idx = np.random.choice(data_domain.shape[0], num, replace=False)
x_train = data_domain[idx, 0:1]
y_train = data_domain[idx, 1:2]
t_train = data_domain[idx, 2:3]
u_train = data_domain[idx, 3:4]
v_train = data_domain[idx, 4:5]
p_train = data_domain[idx, 5:6]
return [x_train, y_train, t_train, u_train, v_train, p_train]
# Parameters to be identified
C1 = dde.Variable(0.0)
C2 = dde.Variable(0.0)
# Define Navier Stokes Equations (Time-dependent PDEs)
def Navier_Stokes_Equation(x, y):
u = y[:, 0:1]
v = y[:, 1:2]
p = y[:, 2:3]
du_x = dde.grad.jacobian(y, x, i=0, j=0)
du_y = dde.grad.jacobian(y, x, i=0, j=1)
du_t = dde.grad.jacobian(y, x, i=0, j=2)
dv_x = dde.grad.jacobian(y, x, i=1, j=0)
dv_y = dde.grad.jacobian(y, x, i=1, j=1)
dv_t = dde.grad.jacobian(y, x, i=1, j=2)
dp_x = dde.grad.jacobian(y, x, i=2, j=0)
dp_y = dde.grad.jacobian(y, x, i=2, j=1)
du_xx = dde.grad.hessian(y, x, component=0, i=0, j=0)
du_yy = dde.grad.hessian(y, x, component=0, i=1, j=1)
dv_xx = dde.grad.hessian(y, x, component=1, i=0, j=0)
dv_yy = dde.grad.hessian(y, x, component=1, i=1, j=1)
continuity = du_x + dv_y
x_momentum = du_t + C1 * (u * du_x + v * du_y) + dp_x - C2 * (du_xx + du_yy)
y_momentum = dv_t + C1 * (u * dv_x + v * dv_y) + dp_y - C2 * (dv_xx + dv_yy)
return [continuity, x_momentum, y_momentum]
# Define Spatio-temporal domain
# Rectangular
Lx_min, Lx_max = 1.0, 8.0
Ly_min, Ly_max = -2.0, 2.0
# Spatial domain: X × Y = [1, 8] × [−2, 2]
space_domain = dde.geometry.Rectangle([Lx_min, Ly_min], [Lx_max, Ly_max])
# Time domain: T = [0, 7]
time_domain = dde.geometry.TimeDomain(0, 7)
# Spatio-temporal domain
geomtime = dde.geometry.GeometryXTime(space_domain, time_domain)
# Get the training data: num = 7000
[ob_x, ob_y, ob_t, ob_u, ob_v, ob_p] = load_training_data(num=7000)
ob_xyt = np.hstack((ob_x, ob_y, ob_t))
observe_u = dde.icbc.PointSetBC(ob_xyt, ob_u, component=0)
observe_v = dde.icbc.PointSetBC(ob_xyt, ob_v, component=1)
# Training datasets and Loss
data = dde.data.TimePDE(
geomtime,
Navier_Stokes_Equation,
[observe_u, observe_v],
num_domain=700,
num_boundary=200,
num_initial=100,
anchors=ob_xyt,
)
# Neural Network setup
layer_size = [3] + [50] * 6 + [3]
activation = "tanh"
initializer = "Glorot uniform"
net = dde.nn.FNN(layer_size, activation, initializer)
model = dde.Model(data, net)
# callbacks for storing results
fnamevar = "variables.dat"
variable = dde.callbacks.VariableValue([C1, C2], period=100, filename=fnamevar)
# Compile, train and save model
model.compile("adam", lr=1e-3, external_trainable_variables=[C1, C2])
loss_history, train_state = model.train(
iterations=10000, callbacks=[variable], display_every=1000, disregard_previous_best=True
)
dde.saveplot(loss_history, train_state, issave=True, isplot=True)
model.compile("adam", lr=1e-4, external_trainable_variables=[C1, C2])
loss_history, train_state = model.train(
iterations=10000, callbacks=[variable], display_every=1000, disregard_previous_best=True
)
dde.saveplot(loss_history, train_state, issave=True, isplot=True)
# model.save(save_path = "./NS_inverse_model/model")
f = model.predict(ob_xyt, operator=Navier_Stokes_Equation)
print("Mean residual:", np.mean(np.absolute(f)))
# Plot Variables:
# reopen saved data using callbacks in fnamevar
lines = open(fnamevar, "r").readlines()
# read output data in fnamevar
Chat = np.array(
[
np.fromstring(
min(re.findall(re.escape("[") + "(.*?)" + re.escape("]"), line), key=len),
sep=",",
)
for line in lines
]
)
l, c = Chat.shape
plt.semilogy(range(0, l * 100, 100), Chat[:, 0], "r-")
plt.semilogy(range(0, l * 100, 100), Chat[:, 1], "k-")
plt.semilogy(range(0, l * 100, 100), np.ones(Chat[:, 0].shape) * C1true, "r--")
plt.semilogy(range(0, l * 100, 100), np.ones(Chat[:, 1].shape) * C2true, "k--")
plt.legend(["C1hat", "C2hat", "True C1", "True C2"], loc="right")
plt.xlabel("Epochs")
plt.title("Variables")
plt.show()
# Plot the velocity distribution of the flow field:
for t in range(0, 8):
[ob_x, ob_y, ob_t, ob_u, ob_v, ob_p] = load_training_data(num=140000)
xyt_pred = np.hstack((ob_x, ob_y, t * np.ones((len(ob_x), 1))))
uvp_pred = model.predict(xyt_pred)
x_pred, y_pred, t_pred = xyt_pred[:, 0], xyt_pred[:, 1], xyt_pred[:, 2]
u_pred, v_pred, p_pred = uvp_pred[:, 0], uvp_pred[:, 1], uvp_pred[:, 2]
x_true = ob_x[ob_t == t]
y_true = ob_y[ob_t == t]
u_true = ob_u[ob_t == t]
fig, ax = plt.subplots(2, 1)
cntr0 = ax[0].tricontourf(x_pred, y_pred, u_pred, levels=80, cmap="rainbow")
cb0 = plt.colorbar(cntr0, ax=ax[0])
cntr1 = ax[1].tricontourf(x_true, y_true, u_true, levels=80, cmap="rainbow")
cb1 = plt.colorbar(cntr1, ax=ax[1])
ax[0].set_title("u-PINN " + "(t=" + str(t) + ")", fontsize=9.5)
ax[0].axis("scaled")
ax[0].set_xlabel("X", fontsize=7.5, family="Arial")
ax[0].set_ylabel("Y", fontsize=7.5, family="Arial")
ax[1].set_title("u-Reference solution " + "(t=" + str(t) + ")", fontsize=9.5)
ax[1].axis("scaled")
ax[1].set_xlabel("X", fontsize=7.5, family="Arial")
ax[1].set_ylabel("Y", fontsize=7.5, family="Arial")
fig.tight_layout()
plt.show()