-
Notifications
You must be signed in to change notification settings - Fork 3
Kuramoto Sivashinky 1D PDE #9
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
base: master
Are you sure you want to change the base?
Changes from 7 commits
b849181
c7beadc
e2f3c6a
dd2f38a
09e0966
7e63fe2
6445c3c
91d439c
9a8a63c
3eb99cd
c1bf87e
1b53abc
33b0987
09f0c92
8c0e7a3
590b7dd
f9af33c
File filter
Filter by extension
Conversations
Jump to
Diff view
Diff view
There are no files selected for viewing
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,49 @@ | ||
% 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. | ||
% | ||
|
||
classdef Canonical < otp.kuramotosivashinsky.KuramotoSivashinskyProblem | ||
|
||
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.L = L; | ||
|
||
h = L/N; | ||
|
||
% exclude the left boundary point as it is identical to the | ||
% right boundary point | ||
x = linspace(h, L, N).'; | ||
|
||
u0 = cos(x/16).*(1+sin(x/16)); | ||
|
||
u0hat = fft(u0); | ||
|
||
tspan = [0, 150]; | ||
|
||
obj = [email protected](tspan, u0hat, params); | ||
|
||
end | ||
end | ||
end |
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,59 @@ | ||
classdef KuramotoSivashinskyProblem < otp.Problem | ||
|
||
methods | ||
function obj = KuramotoSivashinskyProblem(timeSpan, y0, parameters) | ||
[email protected]('Kuramoto-Sivashinsky', [], timeSpan, y0, parameters); | ||
end | ||
end | ||
|
||
methods | ||
function soly = solution2real(soly) | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. To be consistent with Pendulum problem, let's rename to something like |
||
|
||
if isstruct(soly) | ||
soly.y = real(ifft(soly.y)); | ||
else | ||
soly = real(ifft(soly.')).'; | ||
end | ||
|
||
end | ||
|
||
end | ||
|
||
methods (Access = protected) | ||
function onSettingsChanged(obj) | ||
|
||
L = obj.Parameters.L; | ||
|
||
N = obj.NumVars; | ||
|
||
div = L/(2*pi); | ||
|
||
k = (1i*[0:(N/2 - 1), 0, (-N/2 + 1):-1].'/div); | ||
k2 = k.^2; | ||
k4 = k.^4; | ||
Steven-Roberts marked this conversation as resolved.
Show resolved
Hide resolved
|
||
|
||
obj.Rhs = otp.Rhs(@(t, u) otp.kuramotosivashinsky.f(t, u, k, k2, k4), ... | ||
otp.Rhs.FieldNames.Jacobian, ... | ||
@(t, u) otp.kuramotosivashinsky.jac(t,u, k, k2, k4), ... | ||
otp.Rhs.FieldNames.JacobianVectorProduct, ... | ||
@(t, u, v) otp.kuramotosivashinsky.jvp(t, u, v, k, k2, k4), ... | ||
otp.Rhs.FieldNames.JacobianAdjointVectorProduct, ... | ||
@(t, u, v) otp.kuramotosivashinsky.javp(t, u, v, k, k2, k4)); | ||
|
||
end | ||
|
||
function validateNewState(obj, newTimeSpan, newY0, newParameters) | ||
|
||
[email protected](obj, ... | ||
newTimeSpan, newY0, newParameters) | ||
|
||
if mod(numel(newY0), 2) ~= 0 | ||
error('The problem size has to be an even integer.'); | ||
end | ||
|
||
otp.utils.StructParser(newParameters) ... | ||
.checkField('L', 'scalar', 'finite', 'positive'); | ||
|
||
end | ||
end | ||
end |
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,7 @@ | ||
function ut = f(~, u, k, k2, k4) | ||
|
||
u2 = fft(real(ifft(u)).^2); | ||
|
||
ut = - k2.*u - k4.*u -(k/2).*u2; | ||
|
||
end |
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,5 @@ | ||
function j = jac(~, u, k, k2, k4) | ||
|
||
j = -diag(k2) - diag(k4) - diag(k)*ifft(fft(diag(ifft(u))).').'; | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more.
|
||
|
||
end |
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,5 @@ | ||
function jv = javp(~, u, v, k, k2, k4) | ||
|
||
jv = -k2.*v - k4.*v - fft(real(ifft(u)).*ifft(conj(k).*v)); | ||
|
||
end |
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,5 @@ | ||
function jv = jvp(~, u, v, k, k2, k4) | ||
|
||
jv = -k2.*v - k4.*v - k.*fft(real(ifft(u)).*ifft(v)); | ||
|
||
end |
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Replace hardcoded 16 so i.c. is periodic for different L. It might be easier to have x in the range 0 to 2 pi.
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
If I have x in 0 to 2 pi, then I would have to calculate h differently. I will change u0, to be a function of L and use cospi and sinpi. I think that is simpler.