Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Kuramoto Sivashinky 1D PDE #9

Open
wants to merge 17 commits into
base: master
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
45 changes: 45 additions & 0 deletions src/+otp/+kuramotosivashinsky/+presets/Canonical.m
Original file line number Diff line number Diff line change
@@ -0,0 +1,45 @@
classdef Canonical < otp.kuramotosivashinsky.KuramotoSivashinskyProblem
% The Kuramoto-Sivashinky equation is a chaotic problem.
%
% In this particular discretization, we are applying a spectral
% method, therefore the boundary conditions will be chosen to be as cyclic
% on the domain [0, L]. Note that this is different from another typical
% domain of [-L, L]. The larger the L, the more interesting the problem is
% but the more points are required to do a good discretization. The current
% canonical implementation with the size, L, and is used in
%
% Kassam, Aly-Khan, and Lloyd N. Trefethen.
% "Fourth-order time-stepping for stiff PDEs."
% SIAM Journal on Scientific Computing 26, no. 4 (2005): 1214-1233.
%

methods
function obj = Canonical(varargin)

p = inputParser;
addParameter(p, 'Size', 128);
addParameter(p, 'L', 32*pi);
parse(p, varargin{:});

s = p.Results;

N = s.Size;
L = s.L;
params = otp.kuramotosivashinsky.KuramotoSivashinskyParameters;
params.L = L;

% exclude the left boundary point as it is identical to the
% right boundary point
x = linspace(L / N, L, N).';

u0 = cospi(2 * x / L) .* (1 + sinpi(2 * x / L));

u0hat = fft(u0);

tspan = [0, 150];

obj = [email protected](tspan, u0hat, params);

end
end
end
38 changes: 38 additions & 0 deletions src/+otp/+kuramotosivashinsky/+presets/HairerWanner.m
Original file line number Diff line number Diff line change
@@ -0,0 +1,38 @@
classdef HairerWanner < otp.kuramotosivashinsky.KuramotoSivashinskyProblem

methods
function obj = HairerWanner(varargin)

p = inputParser;
addParameter(p, 'Size', 1024);
addParameter(p, 'L', 80 * pi);

parse(p, varargin{:});

s = p.Results;

N = s.Size;
L = s.L;

params.L = L;

% exclude the left boundary point as it is identical to the
% right boundary point
x = linspace(L/N, L, N).';

eta1 = min(x/L, 0.1 - x/L);
eta2 = 20*(x/L - 0.2).*(0.3 - x/L);
eta3 = min(x/L - 0.6, 0.7 - x/L);
eta4 = min(x/L - 0.9, 1 - x/L);

u0 = 16*max(0, max([eta1, eta2, eta3, eta4], [], 2));

u0hat = fft(u0);

tspan = [0, 100];

obj = [email protected](tspan, u0hat, params);

end
end
end
7 changes: 7 additions & 0 deletions src/+otp/+kuramotosivashinsky/KuramotoSivashinskyParameters.m
Original file line number Diff line number Diff line change
@@ -0,0 +1,7 @@
classdef KuramotoSivashinskyParameters
% Parameters for the Kuramoto Sivashinsky problem
properties
% Represents the length-scale of the problem
L
end
end
32 changes: 32 additions & 0 deletions src/+otp/+kuramotosivashinsky/KuramotoSivashinskyProblem.m
Original file line number Diff line number Diff line change
@@ -0,0 +1,32 @@
classdef KuramotoSivashinskyProblem < otp.Problem

methods
function obj = KuramotoSivashinskyProblem(timeSpan, y0, parameters)
[email protected]('Kuramoto-Sivashinsky', [], timeSpan, y0, parameters);
end
end

methods
function soly = convert2grid(~, soly)
soly = abs(ifft(soly));
end
end

methods (Access = protected)
function onSettingsChanged(obj)
N = obj.NumVars;

k = 2 * pi * [0:(N/2 - 1), 0, (-N/2 + 1):-1].' / obj.Parameters.L;
ik = 1i * k;
k24 = k.^2 - k.^4;

obj.RHS = otp.RHS(@(t, u) otp.kuramotosivashinsky.f(t, u, ik, k24), ...
'Jacobian', ...
@(t, u) otp.kuramotosivashinsky.jac(t,u, ik, k24), ...
'JacobianVectorProduct', ...
@(t, u, v) otp.kuramotosivashinsky.jvp(t, u, v, ik, k24), ...
'JacobianAdjointVectorProduct', ...
@(t, u, v) otp.kuramotosivashinsky.javp(t, u, v, ik, k24));
end
end
end
5 changes: 5 additions & 0 deletions src/+otp/+kuramotosivashinsky/f.m
Original file line number Diff line number Diff line change
@@ -0,0 +1,5 @@
function ut = f(~, u, ik, k24)

ut = k24 .* u - (ik/2) .* fft(ifft(u).^2);

end
6 changes: 6 additions & 0 deletions src/+otp/+kuramotosivashinsky/jac.m
Original file line number Diff line number Diff line change
@@ -0,0 +1,6 @@
function j = jac(~, u, ik, k24)

circulantU = toeplitz(u, [u(1); flipud(u(2:end))]);
j = diag(k24) - ik / length(u) .* circulantU;

end
5 changes: 5 additions & 0 deletions src/+otp/+kuramotosivashinsky/javp.m
Original file line number Diff line number Diff line change
@@ -0,0 +1,5 @@
function jv = javp(~, u, v, ik, k24)

jv = k24 .* v + fft(conj(ifft(u)) .* ifft(ik .* v));

end
5 changes: 5 additions & 0 deletions src/+otp/+kuramotosivashinsky/jvp.m
Original file line number Diff line number Diff line change
@@ -0,0 +1,5 @@
function jv = jvp(~, u, v, ik, k24)

jv = k24 .* v - ik .* fft(ifft(u) .* ifft(v));

end