-
Notifications
You must be signed in to change notification settings - Fork 20
/
taylor_green_vortex_2d.m
374 lines (298 loc) · 13.5 KB
/
taylor_green_vortex_2d.m
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
% ----------------------------------------------------------------------- %
% __ __ __ _ __ __ %
% |\/| _ |_ | _ |_ |__| / |_ | \ _ (_ |__) |_ %
% | | (_| |_ | (_| |_) | \__ | |__/ (_) | | \ | %
% %
% ----------------------------------------------------------------------- %
% %
% Author: Alberto Cuoci <[email protected]> %
% CRECK Modeling Group <http://creckmodeling.chem.polimi.it> %
% Department of Chemistry, Materials and Chemical Engineering %
% Politecnico di Milano %
% P.zza Leonardo da Vinci 32, 20133 Milano %
% %
% ----------------------------------------------------------------------- %
% %
% This file is part of Matlab4CFDofRF framework. %
% %
% License %
% %
% Copyright(C) 2019 Alberto Cuoci %
% Matlab4CFDofRF is free software: you can redistribute it and/or %
% modify it under the terms of the GNU General Public License as %
% published by the Free Software Foundation, either version 3 of the %
% License, or (at your option) any later version. %
% %
% Matlab4CFDofRF is distributed in the hope that it will be useful, %
% but WITHOUT ANY WARRANTY; without even the implied warranty of %
% MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the %
% GNU General Public License for more details. %
% %
% You should have received a copy of the GNU General Public License %
% along with Matlab4CRE. If not, see <http://www.gnu.org/licenses/>. %
% %
%-------------------------------------------------------------------------%
% %
% Code: Taylor-Green vortex in 2D
% The Taylor-Green vortex is an exact closed form solution of 2D, %
% incompressible Navier-Stokes eqs. This 2D decaying vortex defined %
% in the square domain, 0-2pi, serves as a benchmark problem %
% for testing and validation of incompressible Navier-Stokes codes. %
% %
% ----------------------------------------------------------------------- %
close all;
clear variables;
% ----------------------------------------------------------------------- %
% User data
% ----------------------------------------------------------------------- %
% Only even numbers of cells are acceptable
nx=60; % number of (physical) cells along x
ny=nx; % number of (physical) cells along y
L=2*pi; % length [m]
nu=1; % kinematic viscosity [m2/s]
tau=1; % total time of simulation [s]
% Parameters for SOR
max_iterations=20000; % maximum number of iterations
beta=1.5; % SOR coefficient
max_error=1e-7; % error for convergence
% Process the grid
h=L/nx; % grid step (uniform grid) [m]
% Initial Solution
t = 0;
[u,v] = AnalyticalSolutionVelocity(nx,ny, h, nu, t);
[p] = AnalyticalSolutionPressure(nx,ny, h, nu, t);
% Time step
umax = sqrt( max(max(u))^2 + ...
max(max(v))^2 ); % maximum velocity [m/s]
sigma = 0.5; % safety factor for time step (stability)
dt_diff=h^2/4/nu; % time step (diffusion stability) [s]
dt_conv=4*nu/umax^2; % time step (convection stability) [s]
dt=sigma*min(dt_diff, dt_conv); % time step (stability) [s]
nsteps=tau/dt; % number of steps
fprintf('Time step: %f\n', dt);
fprintf(' - Diffusion: %f\n', dt_diff);
fprintf(' - Convection: %f\n', dt_conv);
% Grid construction
x=0:h:L; % grid coordinates (x axis)
y=0:h:L; % grid coordinates (y axis)
[X,Y] = meshgrid(x,y); % MATLAB grid
% Temporary velocity fields
ut=u;
vt=v;
% Coefficient for pressure equation
gamma=zeros(nx+2,ny+2)+1/4;
gamma(2,3:ny)=1/3;gamma(nx+1,3:ny)=1/3;gamma(3:nx,2)=1/3;gamma(3:nx,ny+1)=1/3;
gamma(2,2)=1/2;gamma(2,ny+1)=1/2;gamma(nx+1,2)=1/2;gamma(nx+1,ny+1)=1/2;
% ----------------------------------------------------------------------- %
% Solution over time
% ----------------------------------------------------------------------- %
t=0.0;
for is=1:1:nsteps
% Boundary conditions
u(1:nx+1,1)=u(1:nx+1,ny+1); % south wall
u(1:nx+1,ny+2)=u(1:nx+1,2); % north wall
v(1,1:ny+1)=v(nx+1,1:ny+1); % west wall
v(nx+2,1:ny+1)=v(2,1:ny+1); % east wall
% Advection-diffusion equation (predictor)
[ut, vt] = AdvectionDiffusion2D( ut, vt, u, v, nx, ny, h, dt, nu);
% Pressure equation (Poisson)
[p, iter] = Poisson2D( p, ut, vt, gamma, nx, ny, h, dt, ...
beta, max_iterations, max_error );
% Correction on the velocity
u(2:nx,2:ny+1)=ut(2:nx,2:ny+1)-(dt/h)*(p(3:nx+1,2:ny+1)-p(2:nx,2:ny+1));
v(2:nx+1,2:ny)=vt(2:nx+1,2:ny)-(dt/h)*(p(2:nx+1,3:ny+1)-p(2:nx+1,2:ny));
% Print on the screen
if (mod(is,100)==1)
fprintf('Step: %d - Time: %f - Poisson iterations: %d\n', is, t, iter);
% Reconstruct fields
[uu, vv, pp] = FieldReconstruction(u,v,p, nx,ny);
% Surface map: u-velocity
subplot(221);
surface(X,Y,uu');
axis('square'); title('u'); xlabel('x'); ylabel('y');
shading interp; colorbar;
% Surface map: v-velocity
subplot(222);
surface(X,Y,vv');
axis('square'); title('v'); xlabel('x'); ylabel('y');
shading interp; colorbar;
% Surface map: pressure
subplot(223);
surface(X,Y,pp');
axis('square'); title('p'); xlabel('x'); ylabel('y');
shading interp; colorbar;
% Plot: velocity components along the horizontal middle axis
subplot(224);
plot(x,uu(:, round(ny/2)+1));
hold on;
plot(x,vv(:, round(ny/2)+1));
axis('square');
xlabel('x [m]');ylabel('velocity [m/s]'); title('velocity along x-axis');
legend('x-velocity', 'y-velocity');
hold off;
pause(0.0001);
end
% Advance in time
t=t+dt;
end
% ----------------------------------------------------------------------- %
% Final post-processing %
% ----------------------------------------------------------------------- %
% Field reconstruction
[uu, vv, pp] = FieldReconstruction(u,v,p, nx,ny);
% Analytical Solution
[ua,va] = AnalyticalSolutionVelocity(nx,ny,h, nu, t);
[pa] = AnalyticalSolutionPressure(nx,ny,h,nu,t);
[uua, vva, ppa] = FieldReconstruction(ua,va,pa, nx,ny);
% Surface map: u-velocity
subplot(241);
surface(X,Y,uu'); shading interp; colorbar;
axis('square'); title('u (numerical)'); xlabel('x'); ylabel('y');
% Surface map: v-velocity
subplot(242);
surface(X,Y,vv'); shading interp; colorbar;
axis('square'); title('v (numerical)'); xlabel('x'); ylabel('y');
% Surface map: pressure
subplot(243);
surface(X,Y,pp'); shading interp; colorbar;
axis('square'); title('p (numerical)'); xlabel('x'); ylabel('y');
% Surface map: difference pressure
subplot(244);
surface(X,Y,(ppa-pp)'); shading interp; colorbar;
axis('square'); title('delta(p)'); xlabel('x'); ylabel('y');
% Surface map: u-velocity (analytical)
subplot(245);
surface(X,Y,uua'); shading interp; colorbar;
axis('square'); title('u (analytical)'); xlabel('x'); ylabel('y');
% Surface map: v-velocity (analytical)
subplot(246);
surface(X,Y,vva'); shading interp; colorbar;
axis('square'); title('v (analytical)'); xlabel('x'); ylabel('y');
% Surface map: pressure (analytical)
subplot(247);
surface(X,Y,ppa'); shading interp; colorbar;
axis('square'); title('p (analytical)'); xlabel('x'); ylabel('y');
% Streamlines
subplot(248);
sx = [0.25*L 0.75*L 0.25*L 0.75*L];
sy = [0.25*L 0.75*L 0.75*L 0.25*L];
streamline(X,Y,uu',vv',sx,sy)
axis([0 L 0 L], 'square');
title('streamlines'); xlabel('x'); ylabel('y');
% --------------------------------------------------------------------------------------
% Poisson equation solver
% --------------------------------------------------------------------------------------
function [p, iter] = Poisson2D( p, ut, vt, gamma, nx, ny, h, dt, ...
beta, max_iterations, max_error)
for iter=1:max_iterations
for i=2:nx+1
for j=2:ny+1
delta = p(i+1,j)+p(i-1,j)+p(i,j+1)+p(i,j-1);
S = (h/dt)*(ut(i,j)-ut(i-1,j)+vt(i,j)-vt(i,j-1));
p(i,j)=beta*gamma(i,j)*( delta-S )+(1-beta)*p(i,j);
end
end
% Estimate the error
epsilon=0.0;
for i=2:nx+1
for j=2:ny+1
delta = p(i+1,j)+p(i-1,j)+p(i,j+1)+p(i,j-1);
S = (h/dt)*(ut(i,j)-ut(i-1,j)+vt(i,j)-vt(i,j-1));
epsilon=epsilon+abs( p(i,j) - gamma(i,j)*( delta-S ) );
end
end
epsilon = epsilon / (nx*ny);
% Check the error
if (epsilon <= max_error) % stop if converged
break;
end
end
end
% --------------------------------------------------------------------------------------
% Advection-diffusion equation
% --------------------------------------------------------------------------------------
function [ut, vt] = AdvectionDiffusion2D( ut, vt, u, v, nx, ny, h, dt, nu)
% Temporary u-velocity
for i=2:nx
for j=2:ny+1
ue2 = 0.25*( u(i+1,j)+u(i,j) )^2;
uw2 = 0.25*( u(i,j)+u(i-1,j) )^2;
unv = 0.25*( u(i,j+1)+u(i,j) )*( v(i+1,j)+v(i,j) );
usv = 0.25*( u(i,j)+u(i,j-1) )*( v(i+1,j-1)+v(i,j-1) );
A = (ue2-uw2+unv-usv)/h;
D = (nu/h^2)*(u(i+1,j)+u(i-1,j)+u(i,j+1)+u(i,j-1)-4*u(i,j));
ut(i,j)=u(i,j)+dt*(-A+D);
end
end
% Temporary v-velocity
for i=2:nx+1
for j=2:ny
vn2 = 0.25*( v(i,j+1)+v(i,j) )^2;
vs2 = 0.25*( v(i,j)+v(i,j-1) )^2;
veu = 0.25*( u(i,j+1)+u(i,j) )*( v(i+1,j)+v(i,j) );
vwu = 0.25*( u(i-1,j+1)+u(i-1,j) )*( v(i,j)+v(i-1,j) );
A = (vn2 - vs2 + veu - vwu)/h;
D = (nu/h^2)*(v(i+1,j)+v(i-1,j)+v(i,j+1)+v(i,j-1)-4*v(i,j));
vt(i,j)=v(i,j)+dt*(-A+D);
end
end
end
% --------------------------------------------------------------------------------------
% Field reconstruction (for graphical purposes only)
% --------------------------------------------------------------------------------------
function [uu, vv, pp] = FieldReconstruction(u,v,p, nx,ny)
% u-velocity
uu=zeros(nx+1,ny+1);
uu(1:nx+1,1:ny+1)=0.50*(u(1:nx+1,2:ny+2)+u(1:nx+1,1:ny+1));
% v-velocity
vv=zeros(nx+1,ny+1);
vv(1:nx+1,1:ny+1)=0.50*(v(2:nx+2,1:ny+1)+v(1:nx+1,1:ny+1));
% Pressure
pp=zeros(nx+1,ny+1);
p(1:nx+2,1) = p(1:nx+2,ny+1); % south wall
p(1:nx+2,ny+2) = p(1:nx+2,2); % north wall
p(1,1:ny+2) = p(nx+1,1:ny+2); % west wall
p(nx+2, 1:ny+2) = p(2,1:ny+2); % east wall
pp(1:nx+1,1:ny+1)=0.25*(p(1:nx+1,1:ny+1)+p(1:nx+1,2:ny+2)+...
p(2:nx+2,1:ny+1)+p(2:nx+2,2:ny+2) );
end
% --------------------------------------------------------------------------------------
% Analytical solution (velocity)
% --------------------------------------------------------------------------------------
function [u,v] = AnalyticalSolutionVelocity(nx,ny, h, nu, t)
u=zeros(nx+1,ny+2);
v=zeros(nx+2,ny+1);
for i=1:nx+1
for j=1:ny+2
x = (i-1)*h;
y = h/2+(j-2)*h;
u(i,j) = -sin(x)*cos(y)*F(t,nu);
end
end
for i=1:nx+2
for j=1:ny+1
x = h/2+(i-2)*h;
y = (j-1)*h;
v(i,j) = cos(x)*sin(y)*F(t,nu);
end
end
end
% --------------------------------------------------------------------------------------
% Analytical solution (pressure)
% --------------------------------------------------------------------------------------
function [p] = AnalyticalSolutionPressure(nx,ny, h, nu, t)
p=zeros(nx+2,ny+2);
for i=2:nx+1
for j=2:ny+1
x = h/2+(i-2)*h;
y = h/2+(j-2)*h;
p(i,j) = 0.25*(cos(2*x)+cos(2*y))*F(t,nu)^2;
end
end
end
% --------------------------------------------------------------------------------------
% Decaying function (for analytical solution)
% --------------------------------------------------------------------------------------
function f = F(t,nu)
f = exp(-2.*nu*t);
end