-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathse_cfl_condition.m
173 lines (153 loc) · 5.62 KB
/
se_cfl_condition.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
function max_courant = se_cfl_condition(npts, temporal_scheme, nuorder, nuvalue, nutype, numethod)
% se_cfl_condition - Determine the CFL condition (maximum stable Courant
% number) for the spectral element discretization of the advection
% equation with diffusion.
%
% Syntax: se_cfl_condition(npts, temporal_scheme, nuorder, nuvalue, nutype, numethod)
%
% Inputs:
% npts - Order of the spectral element method
% temporal_scheme - Temporal integration scheme from time_discrete_M
% nuorder - Order of diffusion to apply (must be even)
% nuvalue - Non-dimensional diffusion coefficient (scalar or vector)
% nutype - 1 [default] (use laplacian operator for diffusion)
% - 2 (use derivative operator for diffusion)
% - 3 (use Legendre filter applied to shortest wavelength mode)
% numethod - 'inline' [default] (apply diffusion at each RK substage)
% - 'split' (apply diffusion after RK cycle)
%
% Outputs:
% max_courant - Maximum Courant number permitted by the CFL condition.
%
% Example:
% >> se_cfl_condition(4, 'RK3', 4, 0.001);
%
%
% Author: Paul Ullrich
% University of California, Davis
% Email address: [email protected]
% Last revision: 08-Oct-2017
%------------- BEGIN CODE --------------
% Check arguments
if (~isscalar(npts))
error('npts argument must be a scalar');
end
if (~exist('nuvalue'))
if (exist('nuorder'))
error('Diffusion requires both nuorder and nuvalue to be specified');
end
nuorder = 0;
nuvalue = 0;
end
if (mod(nuorder, 2) ~= 0)
error('nuorder argument must be even');
end
if (~isvector(nuvalue))
error('nuvalue argument must be a vector');
end
if (sum(nuvalue < 0) > 0)
disp('WARNING: anti-diffusion being applied in calculation');
end
if (~exist('nutype'))
nutype = 2;
end
if ((nutype ~= 1) && (nutype ~= 2) && (nutype ~= 3))
error('nutype argument must take value 1, 2, or 3');
end
if (~exist('numethod'))
numethod = 'inline';
end
if (strcmp(numethod, 'inline') && strcmp(numethod, 'split'))
error('numethod must be one of inline or split');
end
% Initialize maximum Courant number
max_courant = zeros(size(nuvalue));
% theta is chosen based on search precision
theta = [0:0.001:1] * pi;
% Loop over all nuvalues
for n = 1:length(nuvalue)
% Bisection search to find stable courant number
search_range = [0 3];
while(1)
% Current CFL number to test
dt = 0.5 * sum(search_range);
stable = 1;
% Inline diffusion
if (strcmp(numethod, 'inline'))
if ((nutype == 1) || (nutype == 2))
M = se_advection_matrix_exact(npts, theta, nuorder, nuvalue(n), nutype);
for j = 1:length(theta)
M(:,:,j) = time_discrete_M(temporal_scheme, M(:,:,j), dt);
end
elseif (nutype == 3)
legcoeffs = ones(npts,1);
legcoeffs(npts) = legcoeffs(npts) - nuvalue;
[B1,B2,B3] = se_legfilter_matrices(npts, legcoeffs);
M = se_advection_matrix_exact(npts, theta);
for j = 1:length(theta)
F(:,:,j) = B1 * exp(-1i * theta(j)) + B2 + B3 * exp(1i * theta(j));
M(:,:,j) = (F(:,:,j) - eye(npts-1)) / dt + M(:,:,j);
M(:,:,j) = time_discrete_M(temporal_scheme, M(:,:,j), dt);
end
end
end
% Time split
if (strcmp(numethod, 'split'))
M = se_advection_matrix_exact(npts, theta);
for j = 1:length(theta)
M(:,:,j) = time_discrete_M(temporal_scheme, M(:,:,j), dt);
end
% Diffusion via derivative operator
if (nutype == 1)
[M1,M2,M3] = se_derivative_matrices(npts);
for j = 1:length(theta)
MM = M1 * exp(-1i * theta(j)) + M2 + M3 * exp(1i * theta(j));
M(:,:,j) = (eye(npts-1) - dt * (-1)^(nuorder/2) * nuvalue(n) * MM^(nuorder)) * M(:,:,j);
end
end
% Diffusion via Laplacian operator
if (nutype == 2)
[D1,D2,D3] = se_laplacian_matrices(npts);
for j = 1:length(theta)
D = D1 * exp(-1i * theta(j)) + D2 + D3 * exp(1i * theta(j));
M(:,:,j) = (eye(npts-1) - dt * (-1)^(nuorder/2) * nuvalue(n) * D^(nuorder/2)) * M(:,:,j);
end
end
% Diffusion via Boyd-Vandeven filter
if (nutype == 3)
lfcoeffs = ones(npts,1);
lfcoeffs(npts) = lfcoeffs(npts) - nuvalue;
[B1,B2,B3] = se_legfilter_matrices(npts, lfcoeffs);
for j = 1:length(theta)
F(:,:,j) = B1 * exp(-1i * theta(j)) + B2 + B3 * exp(1i * theta(j));
M(:,:,j) = F(:,:,j) * M(:,:,j);
end
end
end
% Check stability across all wavenumbers
for j = 1:length(theta)
% Compute eigenvalues
evals = eig(M(:,:,j));
for k = 1:length(evals)
if (abs(evals(k)) > 1.0 + 1.0e-14)
stable = 0;
break;
end
end
if (stable == 0)
break;
end
end
if (stable == 1)
search_range(1) = dt;
else
search_range(2) = dt;
end
if (search_range(2) - search_range(1) < 0.00001)
break;
end
end
% Rescale maximum Courant number to the average distance between
% degrees of freedom.
max_courant(n) = (npts-1) * 0.5 * sum(search_range);
end